Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel2.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27
30#include "G4SystemOfUnits.hh"
32#include "G4LossTableManager.hh"
35#include "G4EmParameters.hh"
36#include "G4DNABornAngle.hh"
37#include "G4DeltaAngle.hh"
38#include "G4Exp.hh"
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
42using namespace std;
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45
47 const G4String& nam) :
48G4VEmModel(nam)
49{
50 verboseLevel = 0;
51 // Verbosity scale:
52 // 0 = nothing
53 // 1 = warning for energy non-conservation
54 // 2 = details of energy budget
55 // 3 = calculation of cross sections, file openings, sampling of atoms
56 // 4 = entering in methods
57
58 if (verboseLevel > 0)
59 {
60 G4cout << "Born ionisation model is constructed " << G4endl;
61 }
62
63 // Mark this model as "applicable" for atomic deexcitation
64
66 fAtomDeexcitation = nullptr;
68 fpMolWaterDensity = nullptr;
69 fTableData = nullptr;
70 fLowEnergyLimit = 0;
71 fHighEnergyLimit = 0;
72 fParticleDef = nullptr;
73
74 // Define default angular generator
75
77
78 // Selection of computation method
79
80 fasterCode = false;
81
82 // Selection of stationary mode
83
84 statCode = false;
85
86 // Selection of SP scaling
87
88 spScaling = true;
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
94{
95 delete fTableData;
96 fVecm.clear();
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
102 const G4DataVector& /*cuts*/)
103{
104
105 if (verboseLevel > 3)
106 {
107 G4cout << "Calling G4DNABornIonisationModel2::Initialise()" << G4endl;
108 }
109
110 if(fParticleDef != nullptr && particle != fParticleDef)
111 {
112 G4ExceptionDescription description;
113 description << "You are trying to initialized G4DNABornIonisationModel2 "
114 "for particle "
115 << particle->GetParticleName()
116 << G4endl;
117 description << "G4DNABornIonisationModel2 was already initialised "
118 "for particle:" << fParticleDef->GetParticleName() << G4endl;
119 G4Exception("G4DNABornIonisationModel2::Initialise","bornIonInit",
120 FatalException,description);
121 }
122
123 fParticleDef = particle;
124
125 // Energy limits
126 const char* path = G4FindDataDir("G4LEDATA");
127
128 // ***
129
130 G4String particleName = particle->GetParticleName();
131 std::ostringstream fullFileName;
132 fullFileName << path;
133
134 if(particleName == "e-")
135 {
136 fTableFile = "dna/sigma_ionisation_e_born";
137 fLowEnergyLimit = 11.*eV;
138 fHighEnergyLimit = 1.*MeV;
139
140 if (fasterCode)
141 {
142 fullFileName << "/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat";
143 }
144 else
145 {
146 fullFileName << "/dna/sigmadiff_ionisation_e_born.dat";
147 }
148 }
149 else if(particleName == "proton")
150 {
151 fTableFile = "dna/sigma_ionisation_p_born";
152 fLowEnergyLimit = 500. * keV;
153 fHighEnergyLimit = 100. * MeV;
154
155 if (fasterCode)
156 {
157 fullFileName << "/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat";
158 }
159 else
160 {
161 fullFileName << "/dna/sigmadiff_ionisation_p_born.dat";
162 }
163 }
164
165 // Cross section
166
167 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
168 fTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
169 fTableData->LoadData(fTableFile);
170
171 // Final state
172
173 std::ifstream diffCrossSection(fullFileName.str().c_str());
174
175 if (!diffCrossSection)
176 {
177 G4ExceptionDescription description;
178 description << "Missing data file:" << G4endl << fullFileName.str() << G4endl;
179 G4Exception("G4DNABornIonisationModel2::Initialise","em0003",
180 FatalException,description);
181 }
182
183 // Clear the arrays for re-initialization case (MT mode)
184 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
185
186 fTdummyVec.clear();
187 fVecm.clear();
188
189 for (int j=0; j<5; j++)
190 {
191 fProbaShellMap[j].clear();
192 fDiffCrossSectionData[j].clear();
193 fNrjTransfData[j].clear();
194 }
195
196 //
197
198 fTdummyVec.push_back(0.);
199 while(!diffCrossSection.eof())
200 {
201 G4double tDummy;
202 G4double eDummy;
203 diffCrossSection>>tDummy>>eDummy;
204 if (tDummy != fTdummyVec.back()) fTdummyVec.push_back(tDummy);
205
206 G4double tmp;
207 for (int j=0; j<5; j++)
208 {
209 diffCrossSection>> tmp;
210
211 fDiffCrossSectionData[j][tDummy][eDummy] = tmp;
212
213 if (fasterCode)
214 {
215 fNrjTransfData[j][tDummy][fDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
216 fProbaShellMap[j][tDummy].push_back(fDiffCrossSectionData[j][tDummy][eDummy]);
217 }
218
219 // SI - only if eof is not reached
220 if (!diffCrossSection.eof() && !fasterCode) fDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
221
222 if (!fasterCode) fVecm[tDummy].push_back(eDummy);
223
224 }
225 }
226
227 //
228 SetLowEnergyLimit(fLowEnergyLimit);
229 SetHighEnergyLimit(fHighEnergyLimit);
230
231 if( verboseLevel>0 )
232 {
233 G4cout << "Born ionisation model is initialized " << G4endl
234 << "Energy range: "
235 << LowEnergyLimit() / eV << " eV - "
236 << HighEnergyLimit() / keV << " keV for "
237 << particle->GetParticleName()
238 << G4endl;
239 }
240
241 // Initialize water density pointer
242
243 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
244 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
245
246 if (isInitialised)
247 { return;}
250
251 if (!statCode)
252 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
253
254 isInitialised = true;
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
258
260 const G4ParticleDefinition* particleDefinition,
261 G4double ekin,
262 G4double,
263 G4double)
264{
265 if (verboseLevel > 3)
266 {
267 G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel2"
268 << G4endl;
269 }
270
271 if (particleDefinition != fParticleDef) return 0;
272
273 // Calculate total cross section for model
274
275 G4double sigma=0;
276
277 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
278
279 if (ekin >= fLowEnergyLimit && ekin <= fHighEnergyLimit)
280 {
281 sigma = fTableData->FindValue(ekin);
282
283 // ICRU49 electronic SP scaling - ZF, SI
284
285 if (particleDefinition == G4Proton::ProtonDefinition() && ekin < 70*MeV && spScaling)
286 {
287 G4double A = 1.39241700556072800000E-009 ;
288 G4double B = -8.52610412942622630000E-002 ;
289 sigma = sigma * G4Exp(A*(ekin/eV)+B);
290 }
291 //
292 }
293
294 if (verboseLevel > 2)
295 {
296 G4cout << "__________________________________" << G4endl;
297 G4cout << "G4DNABornIonisationModel2 - XS INFO START" << G4endl;
298 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
299 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
300 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
301 G4cout << "G4DNABornIonisationModel2 - XS INFO END" << G4endl;
302 }
303
304 return sigma*waterDensity;
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
308
309void G4DNABornIonisationModel2::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
310 const G4MaterialCutsCouple* couple,
311 const G4DynamicParticle* particle,
312 G4double,
313 G4double)
314{
315
316 if (verboseLevel > 3)
317 {
318 G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel2"
319 << G4endl;
320 }
321
322 G4double k = particle->GetKineticEnergy();
323
324 if (k >= fLowEnergyLimit && k <= fHighEnergyLimit)
325 {
326 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
327 G4double particleMass = particle->GetDefinition()->GetPDGMass();
328 G4double totalEnergy = k + particleMass;
329 G4double pSquare = k * (totalEnergy + particleMass);
330 G4double totalMomentum = std::sqrt(pSquare);
331
332 G4int ionizationShell = 0;
333
334 if (!fasterCode) ionizationShell = RandomSelect(k);
335
336 // SI: The following protection is necessary to avoid infinite loops :
337 // sigmadiff_ionisation_e_born.dat has non zero partial xs at 18 eV for shell 3 (ionizationShell ==2)
338 // sigmadiff_cumulated_ionisation_e_born.dat has zero cumulated partial xs at 18 eV for shell 3 (ionizationShell ==2)
339 // this is due to the fact that the max allowed transfered energy is (18+10.79)/2=17.025 eV and only transfered energies
340 // strictly above this value have non zero partial xs in sigmadiff_ionisation_e_born.dat (starting at trans = 17.12 eV)
341
342 if (fasterCode)
343 do
344 {
345 ionizationShell = RandomSelect(k);
346 } while (k<19*eV && ionizationShell==2 && particle->GetDefinition()==G4Electron::ElectronDefinition());
347
348 G4double secondaryKinetic=-1000*eV;
349
350 if (!fasterCode)
351 {
352 secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
353 }
354 else
355 {
356 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
357 }
358
359 G4int Z = 8;
360
361 G4ThreeVector deltaDirection =
362 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
363 Z, ionizationShell,
364 couple->GetMaterial());
365
366 if (secondaryKinetic>0)
367 {
368 auto dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
369 fvect->push_back(dp);
370 }
371
373 {
374 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
375
376 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
377 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
378 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
379 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
380 finalPx /= finalMomentum;
381 finalPy /= finalMomentum;
382 finalPz /= finalMomentum;
383
384 G4ThreeVector direction;
385 direction.set(finalPx,finalPy,finalPz);
386
387 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit());
388 }
389
390 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
391
392 // AM: sample deexcitation
393 // here we assume that H_{2}O electronic levels are the same as Oxygen.
394 // this can be considered true with a rough 10% error in energy on K-shell,
395
396 std::size_t secNumberInit = 0;
397 std::size_t secNumberFinal = 0;
398
399 G4double bindingEnergy = 0;
400 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
401
402 // SI: additional protection if tcs interpolation method is modified
403 if (k<bindingEnergy) return;
404 //
405
406 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
407
408 // SI: only atomic deexcitation from K shell is considered
409 if((fAtomDeexcitation != nullptr) && ionizationShell == 4)
410 {
411 const G4AtomicShell* shell =
412 fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
413 secNumberInit = fvect->size();
414 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
415 secNumberFinal = fvect->size();
416
417 if(secNumberFinal > secNumberInit)
418 {
419 for (std::size_t i=secNumberInit; i<secNumberFinal; ++i)
420 {
421 //Check if there is enough residual energy
422 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
423 {
424 //Ok, this is a valid secondary: keep it
425 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
426 }
427 else
428 {
429 //Invalid secondary: not enough energy to create it!
430 //Keep its energy in the local deposit
431 delete (*fvect)[i];
432 (*fvect)[i]=nullptr;
433 }
434 }
435 }
436
437 }
438
439 //This should never happen
440 if(bindingEnergy < 0.0)
441 G4Exception("G4DNAEmfietzoglouIonisatioModel1::SampleSecondaries()",
442 "em2050",FatalException,"Negative local energy deposit");
443
444 //bindingEnergy has been decreased
445 //by the amount of energy taken away by deexc. products
446 if (!statCode)
447 {
448 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
449 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
450 }
451 else
452 {
453 fParticleChangeForGamma->SetProposedKineticEnergy(k);
454 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy);
455 }
456
457 // TEST //////////////////////////
458 // if (secondaryKinetic<0) abort();
459 // if (scatteredEnergy<0) abort();
460 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
461 // if (k-scatteredEnergy<0) abort();
462 /////////////////////////////////
463
464 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
466 ionizationShell,
467 theIncomingTrack);
468 }
469}
470
471//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
472
473G4double G4DNABornIonisationModel2::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
474 G4double k,
475 G4int shell)
476{
477 // G4cout << "*** SLOW computation for " << " " << particleDefinition->GetParticleName() << G4endl;
478
479 if (particleDefinition == G4Electron::ElectronDefinition())
480 {
481 G4double maximumEnergyTransfer = 0.;
482 if ((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
483 maximumEnergyTransfer = k;
484 else
485 maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell)) / 2.;
486
487 // SI : original method
488 /*
489 G4double crossSectionMaximum = 0.;
490 for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
491 {
492 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
493 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
494 }
495 */
496
497 // SI : alternative method
498 G4double crossSectionMaximum = 0.;
499
500 G4double minEnergy = waterStructure.IonisationEnergy(shell);
501 G4double maxEnergy = maximumEnergyTransfer;
502 G4int nEnergySteps = 50;
503
504 G4double value(minEnergy);
505 G4double stpEnergy(std::pow(maxEnergy / value,
506 1. / static_cast<G4double>(nEnergySteps - 1)));
507 G4int step(nEnergySteps);
508 while (step > 0)
509 {
510 step--;
511 G4double differentialCrossSection =
512 DifferentialCrossSection(particleDefinition,
513 k / eV,
514 value / eV,
515 shell);
516 if (differentialCrossSection >= crossSectionMaximum)
517 crossSectionMaximum = differentialCrossSection;
518 value *= stpEnergy;
519 }
520 //
521
522 G4double secondaryElectronKineticEnergy = 0.;
523 do
524 {
525 secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
526 } while(G4UniformRand()*crossSectionMaximum >
527 DifferentialCrossSection(particleDefinition, k/eV,
528 (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
529
530 return secondaryElectronKineticEnergy;
531
532 }
533
534 if (particleDefinition == G4Proton::ProtonDefinition())
535 {
536 G4double maximumKineticEnergyTransfer = 4.
537 * (electron_mass_c2 / proton_mass_c2) * k;
538
539 G4double crossSectionMaximum = 0.;
540 for (G4double value = waterStructure.IonisationEnergy(shell);
541 value <= 4. * waterStructure.IonisationEnergy(shell); value += 0.1 * eV)
542 {
543 G4double differentialCrossSection =
544 DifferentialCrossSection(particleDefinition,
545 k / eV,
546 value / eV,
547 shell);
548 if (differentialCrossSection >= crossSectionMaximum)
549 crossSectionMaximum = differentialCrossSection;
550 }
551
552 G4double secondaryElectronKineticEnergy = 0.;
553 do
554 {
555 secondaryElectronKineticEnergy = G4UniformRand()* maximumKineticEnergyTransfer;
556 } while(G4UniformRand()*crossSectionMaximum >=
557 DifferentialCrossSection(particleDefinition, k/eV,
558 (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
559
560 return secondaryElectronKineticEnergy;
561 }
562
563 return 0;
564}
565
566//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
567
568// The following section is not used anymore but is kept for memory
569// GetAngularDistribution()->SampleDirectionForShell is used instead
570
571/*
572 void G4DNABornIonisationModel2::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
573 G4double k,
574 G4double secKinetic,
575 G4double & cosTheta,
576 G4double & phi )
577 {
578 if (particleDefinition == G4Electron::ElectronDefinition())
579 {
580 phi = twopi * G4UniformRand();
581 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
582 else if (secKinetic <= 200.*eV)
583 {
584 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
585 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
586 }
587 else
588 {
589 G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
590 cosTheta = std::sqrt(1.-sin2O);
591 }
592 }
593
594 else if (particleDefinition == G4Proton::ProtonDefinition())
595 {
596 G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
597 phi = twopi * G4UniformRand();
598
599 // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
600
601 // Restriction below 100 eV from Emfietzoglou (2000)
602
603 if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
604 else cosTheta = (2.*G4UniformRand())-1.;
605
606 }
607 }
608 */
609
610//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
612 G4double k,
613 G4double energyTransfer,
614 G4int ionizationLevelIndex)
615{
616 G4double sigma = 0.;
617
618 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)/eV)
619 {
620 G4double valueT1 = 0;
621 G4double valueT2 = 0;
622 G4double valueE21 = 0;
623 G4double valueE22 = 0;
624 G4double valueE12 = 0;
625 G4double valueE11 = 0;
626
627 G4double xs11 = 0;
628 G4double xs12 = 0;
629 G4double xs21 = 0;
630 G4double xs22 = 0;
631
632 // Protection against out of boundary access - proton case : 100 MeV
633 if (k==fTdummyVec.back()) k=k*(1.-1e-12);
634 //
635
636 // k should be in eV and energy transfer eV also
637
638 auto t2 = std::upper_bound(fTdummyVec.begin(),
639 fTdummyVec.end(),
640 k);
641
642 auto t1 = t2 - 1;
643
644 // SI : the following condition avoids situations where energyTransfer >last vector element
645
646 if (energyTransfer <= fVecm[(*t1)].back()
647 && energyTransfer <= fVecm[(*t2)].back())
648 {
649 auto e12 = std::upper_bound(fVecm[(*t1)].begin(),
650 fVecm[(*t1)].end(),
651 energyTransfer);
652 auto e11 = e12 - 1;
653
654 auto e22 = std::upper_bound(fVecm[(*t2)].begin(),
655 fVecm[(*t2)].end(),
656 energyTransfer);
657 auto e21 = e22 - 1;
658
659 valueT1 = *t1;
660 valueT2 = *t2;
661 valueE21 = *e21;
662 valueE22 = *e22;
663 valueE12 = *e12;
664 valueE11 = *e11;
665
666 xs11 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
667 xs12 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
668 xs21 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
669 xs22 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
670
671 }
672
673 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
674 if (xsProduct != 0.)
675 {
676 sigma = QuadInterpolator(valueE11,
677 valueE12,
678 valueE21,
679 valueE22,
680 xs11,
681 xs12,
682 xs21,
683 xs22,
684 valueT1,
685 valueT2,
686 k,
687 energyTransfer);
688 }
689 }
690
691 return sigma;
692}
693
694//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
695
696G4double G4DNABornIonisationModel2::Interpolate(G4double e1,
697 G4double e2,
698 G4double e,
699 G4double xs1,
700 G4double xs2)
701{
702 G4double value = 0.;
703
704 // Log-log interpolation by default
705
706 if (e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
707 && !fasterCode)
708 {
709 G4double a = (std::log10(xs2) - std::log10(xs1))
710 / (std::log10(e2) - std::log10(e1));
711 G4double b = std::log10(xs2) - a * std::log10(e2);
712 G4double sigma = a * std::log10(e) + b;
713 value = (std::pow(10., sigma));
714 }
715
716 // Switch to lin-lin interpolation
717 /*
718 if ((e2-e1)!=0)
719 {
720 G4double d1 = xs1;
721 G4double d2 = xs2;
722 value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
723 }
724 */
725
726 // Switch to log-lin interpolation for faster code
727 if ((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
728 {
729 G4double d1 = std::log10(xs1);
730 G4double d2 = std::log10(xs2);
731 value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
732 }
733
734 // Switch to lin-lin interpolation for faster code
735 // in case one of xs1 or xs2 (=cum proba) value is zero
736
737 if ((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
738 {
739 G4double d1 = xs1;
740 G4double d2 = xs2;
741 value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
742 }
743
744 /*
745 G4cout
746 << e1 << " "
747 << e2 << " "
748 << e << " "
749 << xs1 << " "
750 << xs2 << " "
751 << value
752 << G4endl;
753 */
754
755 return value;
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
759
760G4double G4DNABornIonisationModel2::QuadInterpolator(G4double e11,
761 G4double e12,
762 G4double e21,
763 G4double e22,
764 G4double xs11,
765 G4double xs12,
766 G4double xs21,
767 G4double xs22,
768 G4double t1,
769 G4double t2,
770 G4double t,
771 G4double e)
772{
773 G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
774 G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
775 G4double value = Interpolate(t1,
776 t2,
777 t,
778 interpolatedvalue1,
779 interpolatedvalue2);
780
781 return value;
782}
783
784//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
785
787 G4int level,
788 const G4ParticleDefinition* /*particle*/,
789 G4double kineticEnergy)
790{
791 return fTableData->GetComponent(level)->FindValue(kineticEnergy);
792}
793
794//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
795
796G4int G4DNABornIonisationModel2::RandomSelect(G4double k)
797{
798 G4int level = 0;
799
800 auto valuesBuffer = new G4double[fTableData->NumberOfComponents()];
801 const auto n = (G4int)fTableData->NumberOfComponents();
802 G4int i(n);
803 G4double value = 0.;
804
805 while (i > 0)
806 {
807 i--;
808 valuesBuffer[i] = fTableData->GetComponent(i)->FindValue(k);
809 value += valuesBuffer[i];
810 }
811
812 value *= G4UniformRand();
813
814 i = n;
815
816 while (i > 0)
817 {
818 i--;
819
820 if (valuesBuffer[i] > value)
821 {
822 delete[] valuesBuffer;
823 return i;
824 }
825 value -= valuesBuffer[i];
826 }
827
828
829 delete[] valuesBuffer;
830
831 return level;
832}
833
834//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
835
836G4double G4DNABornIonisationModel2::RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition* particleDefinition,
837 G4double k,
838 G4int shell)
839{
840 // G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
841
842 G4double secondaryElectronKineticEnergy = 0.;
843
844 G4double random = G4UniformRand();
845
846 secondaryElectronKineticEnergy = TransferedEnergy(particleDefinition,
847 k / eV,
848 shell,
849 random) * eV
850 - waterStructure.IonisationEnergy(shell);
851
852 // G4cout << TransferedEnergy(particleDefinition, k/eV, shell, random) << G4endl;
853 if (secondaryElectronKineticEnergy < 0.)
854 return 0.;
855 //
856
857 return secondaryElectronKineticEnergy;
858}
859
860//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
861
863 G4double k,
864 G4int ionizationLevelIndex,
865 G4double random)
866{
867
868 G4double nrj = 0.;
869
870 G4double valueK1 = 0;
871 G4double valueK2 = 0;
872 G4double valuePROB21 = 0;
873 G4double valuePROB22 = 0;
874 G4double valuePROB12 = 0;
875 G4double valuePROB11 = 0;
876
877 G4double nrjTransf11 = 0;
878 G4double nrjTransf12 = 0;
879 G4double nrjTransf21 = 0;
880 G4double nrjTransf22 = 0;
881
882 // Protection against out of boundary access - proton case : 100 MeV
883 if (k==fTdummyVec.back()) k=k*(1.-1e-12);
884 //
885
886 // k should be in eV
887 auto k2 = std::upper_bound(fTdummyVec.begin(),
888 fTdummyVec.end(),
889 k);
890 auto k1 = k2 - 1;
891
892 /*
893 G4cout << "----> k=" << k
894 << " " << *k1
895 << " " << *k2
896 << " " << random
897 << " " << ionizationLevelIndex
898 << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
899 << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
900 << G4endl;
901 */
902
903 // SI : the following condition avoids situations where random >last vector element
904 if (random <= fProbaShellMap[ionizationLevelIndex][(*k1)].back()
905 && random <= fProbaShellMap[ionizationLevelIndex][(*k2)].back())
906 {
907 auto prob12 =
908 std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
909 fProbaShellMap[ionizationLevelIndex][(*k1)].end(),
910 random);
911
912 auto prob11 = prob12 - 1;
913
914 auto prob22 =
915 std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
916 fProbaShellMap[ionizationLevelIndex][(*k2)].end(),
917 random);
918
919 auto prob21 = prob22 - 1;
920
921 valueK1 = *k1;
922 valueK2 = *k2;
923 valuePROB21 = *prob21;
924 valuePROB22 = *prob22;
925 valuePROB12 = *prob12;
926 valuePROB11 = *prob11;
927
928 /*
929 G4cout << " " << random << " " << valuePROB11 << " "
930 << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
931 */
932
933 nrjTransf11 = fNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
934 nrjTransf12 = fNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
935 nrjTransf21 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
936 nrjTransf22 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
937
938 /*
939 G4cout << " " << ionizationLevelIndex << " "
940 << random << " " <<valueK1 << " " << valueK2 << G4endl;
941
942 G4cout << " " << random << " " << nrjTransf11 << " "
943 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
944 */
945
946 }
947 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
948 if (random > fProbaShellMap[ionizationLevelIndex][(*k1)].back())
949 {
950 auto prob22 =
951 std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
952 fProbaShellMap[ionizationLevelIndex][(*k2)].end(),
953 random);
954
955 auto prob21 = prob22 - 1;
956
957 valueK1 = *k1;
958 valueK2 = *k2;
959 valuePROB21 = *prob21;
960 valuePROB22 = *prob22;
961
962 // G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
963
964 nrjTransf21 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
965 nrjTransf22 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
966
967 G4double interpolatedvalue2 = Interpolate(valuePROB21,
968 valuePROB22,
969 random,
970 nrjTransf21,
971 nrjTransf22);
972
973 // zeros are explicitly set
974
975 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
976
977 /*
978 G4cout << " " << ionizationLevelIndex << " "
979 << random << " " <<valueK1 << " " << valueK2 << G4endl;
980
981 G4cout << " " << random << " " << nrjTransf11 << " "
982 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
983
984 G4cout << "ici" << " " << value << G4endl;
985 */
986
987 return value;
988 }
989
990 // End electron and proton cases
991
992 G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21
993 * nrjTransf22;
994 //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
995
996 if (nrjTransfProduct != 0.)
997 {
998 nrj = QuadInterpolator(valuePROB11,
999 valuePROB12,
1000 valuePROB21,
1001 valuePROB22,
1002 nrjTransf11,
1003 nrjTransf12,
1004 nrjTransf21,
1005 nrjTransf22,
1006 valueK1,
1007 valueK2,
1008 k,
1009 random);
1010 }
1011 // G4cout << nrj << endl;
1012
1013 return nrj;
1014}
@ eIonizedMolecule
G4double B(G4double temperature)
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4ThreeVector G4ParticleMomentum
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double) override
G4ParticleChangeForGamma * fParticleChangeForGamma
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4DNABornIonisationModel2(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
size_t NumberOfComponents() const override
const G4VEMDataSet * GetComponent(G4int componentId) const override
static G4DNAMolecularMaterial * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4bool DNAStationary() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
std::size_t GetIndex() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
void SetHighEnergyLimit(G4double)
G4VEmAngularDistribution * GetAngularDistribution()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)