Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel1.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 "G4EmParameters.hh"
33#include "G4LossTableManager.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
65 fAtomDeexcitation = nullptr;
67 fpMolWaterDensity = nullptr;
68
69 // Define default angular generator
71
72 fasterCode = G4EmParameters::Instance()->DNAFast();
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78{
79 // Cross section
80
81 std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
82 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
83 {
84 G4DNACrossSectionDataSet* table = pos->second;
85 delete table;
86 }
87
88 // Final state
89
90 eVecm.clear();
91 pVecm.clear();
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95
97 const G4DataVector& /*cuts*/)
98{
99
100 if (verboseLevel > 3)
101 {
102 G4cout << "Calling G4DNABornIonisationModel1::Initialise()" << G4endl;
103 }
104
105 // Energy limits
106
107 G4String fileElectron("dna/sigma_ionisation_e_born");
108 G4String fileProton("dna/sigma_ionisation_p_born");
109
112
113 G4String electron;
114 G4String proton;
115
116 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
117
118 const char *path = G4FindDataDir("G4LEDATA");
119
120 // *** ELECTRON
121
122 electron = electronDef->GetParticleName();
123
124 tableFile[electron] = fileElectron;
125
126 lowEnergyLimit[electron] = 11. * eV;
127 highEnergyLimit[electron] = 1. * MeV;
128
129 // Cross section
130
131 auto tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
132 tableE->LoadData(fileElectron);
133
134 tableData[electron] = tableE;
135
136 // Final state
137
138 std::ostringstream eFullFileName;
139
140 if (fasterCode) eFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat";
141 if (!fasterCode) eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
142
143 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
144
145 if (!eDiffCrossSection)
146 {
147 if (fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
148 FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat");
149
150 if (!fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
151 FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_born.dat");
152 }
153
154 // Clear the arrays for re-initialization case (MT mode)
155 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
156
157 eTdummyVec.clear();
158 pTdummyVec.clear();
159
160 eVecm.clear();
161 pVecm.clear();
162
163 for (G4int j=0; j<5; j++)
164 {
165 eProbaShellMap[j].clear();
166 pProbaShellMap[j].clear();
167
168 eDiffCrossSectionData[j].clear();
169 pDiffCrossSectionData[j].clear();
170
171 eNrjTransfData[j].clear();
172 pNrjTransfData[j].clear();
173 }
174
175 //
176
177 eTdummyVec.push_back(0.);
178 while(!eDiffCrossSection.eof())
179 {
180 G4double tDummy;
181 G4double eDummy;
182 eDiffCrossSection>>tDummy>>eDummy;
183 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
184
185 G4double tmp;
186 for (G4int j=0; j<5; j++)
187 {
188 eDiffCrossSection>> tmp;
189
190 eDiffCrossSectionData[j][tDummy][eDummy] = tmp;
191
192 if (fasterCode)
193 {
194 eNrjTransfData[j][tDummy][eDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
195 eProbaShellMap[j][tDummy].push_back(eDiffCrossSectionData[j][tDummy][eDummy]);
196 }
197
198 // SI - only if eof is not reached
199 if (!eDiffCrossSection.eof() && !fasterCode) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
200
201 if (!fasterCode) eVecm[tDummy].push_back(eDummy);
202
203 }
204 }
205
206 // *** PROTON
207
208 proton = protonDef->GetParticleName();
209
210 tableFile[proton] = fileProton;
211
212 lowEnergyLimit[proton] = 500. * keV;
213 highEnergyLimit[proton] = 100. * MeV;
214
215 // Cross section
216
217 auto tableP = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
218 tableP->LoadData(fileProton);
219
220 tableData[proton] = tableP;
221
222 // Final state
223
224 std::ostringstream pFullFileName;
225
226 if (fasterCode) pFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat";
227
228 if (!fasterCode) pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
229
230 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
231
232 if (!pDiffCrossSection)
233 {
234 if (fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
235 FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat");
236
237 if (!fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
238 FatalException,"Missing data file:/dna/sigmadiff_ionisation_p_born.dat");
239 }
240
241 pTdummyVec.push_back(0.);
242 while(!pDiffCrossSection.eof())
243 {
244 G4double tDummy;
245 G4double eDummy;
246 pDiffCrossSection>>tDummy>>eDummy;
247 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
248 for (G4int j=0; j<5; j++)
249 {
250 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
251
252 if (fasterCode)
253 {
254 pNrjTransfData[j][tDummy][pDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
255 pProbaShellMap[j][tDummy].push_back(pDiffCrossSectionData[j][tDummy][eDummy]);
256 }
257
258 // SI - only if eof is not reached !
259 if (!pDiffCrossSection.eof() && !fasterCode) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
260
261 if (!fasterCode) pVecm[tDummy].push_back(eDummy);
262 }
263 }
264
265 //
266
267 if (particle==electronDef)
268 {
269 SetLowEnergyLimit(lowEnergyLimit[electron]);
270 SetHighEnergyLimit(highEnergyLimit[electron]);
271 }
272
273 if (particle==protonDef)
274 {
275 SetLowEnergyLimit(lowEnergyLimit[proton]);
276 SetHighEnergyLimit(highEnergyLimit[proton]);
277 }
278
279 if( verboseLevel>0 )
280 {
281 G4cout << "Born ionisation model is initialized " << G4endl
282 << "Energy range: "
283 << LowEnergyLimit() / eV << " eV - "
284 << HighEnergyLimit() / keV << " keV for "
285 << particle->GetParticleName()
286 << G4endl;
287 }
288
289 if (isInitialised) { return; }
291
292 // Initialize water density pointer
293 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
294 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
295
297
298 // AD
299 if (!statCode)
300 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
301
302 // chemistry
304 if (chem->IsChemistryActivated()) {
305 fChemistry = chem;
306 }
307 isInitialised = true;
308}
309
310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311
313{
314 fTrack = track;
315}
316
317//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
318
320 const G4ParticleDefinition* particleDefinition,
321 G4double ekin,
322 G4double,
323 G4double)
324{
325 if (verboseLevel > 3)
326 {
327 G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel1"
328 << G4endl;
329 }
330
331 if (
332 particleDefinition != G4Proton::ProtonDefinition()
333 &&
334 particleDefinition != G4Electron::ElectronDefinition()
335 )
336
337 return 0;
338
339 // Calculate total cross section for model
340
341 G4double lowLim = 0;
342 G4double highLim = 0;
343 G4double sigma=0;
344
345 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
346
347 const G4String& particleName = particleDefinition->GetParticleName();
348
349 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
350 pos1 = lowEnergyLimit.find(particleName);
351 if (pos1 != lowEnergyLimit.end())
352 {
353 lowLim = pos1->second;
354 }
355
356 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
357 pos2 = highEnergyLimit.find(particleName);
358 if (pos2 != highEnergyLimit.end())
359 {
360 highLim = pos2->second;
361 }
362
363 if (ekin >= lowLim && ekin <= highLim)
364 {
365 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
366 pos = tableData.find(particleName);
367
368 if (pos != tableData.end())
369 {
370 G4DNACrossSectionDataSet* table = pos->second;
371 if (table != nullptr)
372 {
373 sigma = table->FindValue(ekin);
374
375 // ICRU49 electronic SP scaling - ZF, SI
376
377 if (particleDefinition == G4Proton::ProtonDefinition() && ekin < 70*MeV && spScaling)
378 {
379 G4double A = 1.39241700556072800000E-009 ;
380 G4double B = -8.52610412942622630000E-002 ;
381 sigma = sigma * G4Exp(A*(ekin/eV)+B);
382 }
383 //
384
385 }
386 }
387 else
388 {
389 G4Exception("G4DNABornIonisationModel1::CrossSectionPerVolume","em0002",
390 FatalException,"Model not applicable to particle type.");
391 }
392 }
393
394 if (verboseLevel > 2)
395 {
396 G4cout << "__________________________________" << G4endl;
397 G4cout << "G4DNABornIonisationModel1 - XS INFO START" << G4endl;
398 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
399 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
400 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
401 G4cout << "G4DNABornIonisationModel1 - XS INFO END" << G4endl;
402 }
403
404 return sigma*waterDensity;
405}
406
407//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
408
409void G4DNABornIonisationModel1::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
410 const G4MaterialCutsCouple* couple,
411 const G4DynamicParticle* particle,
412 G4double,
413 G4double)
414{
415
416 if (verboseLevel > 3)
417 {
418 G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel1"
419 << G4endl;
420 }
421
422 G4double lowLim = 0;
423 G4double highLim = 0;
424
425 G4double k = particle->GetKineticEnergy();
426
427 const G4String& particleName = particle->GetDefinition()->GetParticleName();
428
429 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
430 pos1 = lowEnergyLimit.find(particleName);
431
432 if (pos1 != lowEnergyLimit.end())
433 {
434 lowLim = pos1->second;
435 }
436
437 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
438 pos2 = highEnergyLimit.find(particleName);
439
440 if (pos2 != highEnergyLimit.end())
441 {
442 highLim = pos2->second;
443 }
444
445 if (k >= lowLim && k <= highLim)
446 {
447 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
448 G4double particleMass = particle->GetDefinition()->GetPDGMass();
449 G4double totalEnergy = k + particleMass;
450 G4double pSquare = k * (totalEnergy + particleMass);
451 G4double totalMomentum = std::sqrt(pSquare);
452
453 G4int ionizationShell = 0;
454
455 // SI: The following protection is necessary to avoid infinite loops :
456 // sigmadiff_ionisation_e_born.dat has non zero partial xs at 18 eV for shell 3 (ionizationShell ==2)
457 // sigmadiff_cumulated_ionisation_e_born.dat has zero cumulated partial xs at 18 eV for shell 3 (ionizationShell ==2)
458 // this is due to the fact that the max allowed transfered energy is (18+10.79)/2=17.025 eV and only transfered energies
459 // strictly above this value have non zero partial xs in sigmadiff_ionisation_e_born.dat (starting at trans = 17.12 eV)
460
461 if (!fasterCode) {
462 ionizationShell = RandomSelect(k,particleName);
463 } else {
464 do {
465 ionizationShell = RandomSelect(k,particleName);
466 } while (k<19*eV && ionizationShell==2 && particle->GetDefinition()==G4Electron::ElectronDefinition());
467 }
468 G4double bindingEnergy = 0;
469 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
470
471 // SI: additional protection if tcs interpolation method is modified
472 if (k<bindingEnergy) return;
473 //
474
475 G4double secondaryKinetic=-1000*eV;
476
477 if (!fasterCode)
478 {
479 secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
480 }
481 else
482 {
483 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
484 }
485 G4int Z = 8;
486
487 G4ThreeVector deltaDirection =
488 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
489 Z, ionizationShell,
490 couple->GetMaterial());
491
492 if (secondaryKinetic>0)
493 {
494 auto dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
495 fvect->push_back(dp);
496 }
497
499 {
500 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
501
502 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
503 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
504 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
505 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
506 finalPx /= finalMomentum;
507 finalPy /= finalMomentum;
508 finalPz /= finalMomentum;
509
510 G4ThreeVector direction;
511 direction.set(finalPx,finalPy,finalPz);
512
513 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit());
514 }
515
516 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
517
518 // AM: sample deexcitation
519 // here we assume that H_{2}O electronic levels are the same as Oxygen.
520 // this can be considered true with a rough 10% error in energy on K-shell,
521
522 std::size_t secNumberInit = 0;// need to know at a certain point the energy of secondaries
523 std::size_t secNumberFinal = 0;// So I'll make the diference and then sum the energies
524
525 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
526
527 // SI: only atomic deexcitation from K shell is considered
528 if((fAtomDeexcitation != nullptr) && ionizationShell == 4)
529 {
530 const G4AtomicShell* shell =
531 fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
532 secNumberInit = fvect->size();
533 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
534 secNumberFinal = fvect->size();
535
536 //TEST
537 //G4cout << "ionizationShell=" << ionizationShell<< G4endl;
538 //G4cout << "bindingEnergy=" << bindingEnergy/eV<< G4endl;
539
540 if(secNumberFinal > secNumberInit)
541 {
542 for (std::size_t i=secNumberInit; i<secNumberFinal; ++i)
543 {
544 //Check if there is enough residual energy
545 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
546 {
547 //Ok, this is a valid secondary: keep it
548 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
549 //G4cout << "--deex nrj=" << ((*fvect)[i])->GetKineticEnergy()/eV
550 //<< G4endl;
551 }
552 else
553 {
554 //Invalid secondary: not enough energy to create it!
555 //Keep its energy in the local deposit
556 delete (*fvect)[i];
557 (*fvect)[i]=nullptr;
558 }
559 }
560 }
561
562 //TEST
563 //G4cout << "k=" << k/eV<< G4endl;
564 //G4cout << "secondaryKinetic=" << secondaryKinetic/eV<< G4endl;
565 //G4cout << "scatteredEnergy=" << scatteredEnergy/eV<< G4endl;
566 //G4cout << "deposited energy=" << bindingEnergy/eV<< G4endl;
567 //
568
569 }
570
571 //This should never happen
572 if(bindingEnergy < 0.0)
573 G4Exception("G4DNABornIonisatioModel1::SampleSecondaries()",
574 "em2050",FatalException,"Negative local energy deposit");
575
576 //bindingEnergy has been decreased
577 //by the amount of energy taken away by deexc. products
578 if (!statCode)
579 {
580 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
581 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
582 }
583 else
584 {
585 fParticleChangeForGamma->SetProposedKineticEnergy(k);
586 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy);
587 }
588
589 // create radical
590 if (nullptr != fChemistry) {
591 fChemistry->CreateWaterMolecule(eIonizedMolecule, ionizationShell, fTrack);
592 }
593 }
594}
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
597
598G4double G4DNABornIonisationModel1::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
599 G4double k,
600 G4int shell)
601{
602 // G4cout << "*** SLOW computation for " << " " << particleDefinition->GetParticleName() << G4endl;
603
604 if (particleDefinition == G4Electron::ElectronDefinition())
605 {
606 G4double maximumEnergyTransfer = 0.;
607 if ((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
608 maximumEnergyTransfer = k;
609 else
610 maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell)) / 2.;
611
612 // SI : original method
613 /*
614 G4double crossSectionMaximum = 0.;
615 for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
616 {
617 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
618 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
619 }
620 */
621
622 // SI : alternative method
623 G4double crossSectionMaximum = 0.;
624
625 G4double minEnergy = waterStructure.IonisationEnergy(shell);
626 G4double maxEnergy = maximumEnergyTransfer;
627 G4int nEnergySteps = 50;
628
629 G4double value(minEnergy);
630 G4double stpEnergy(std::pow(maxEnergy / value,
631 1. / static_cast<G4double>(nEnergySteps - 1)));
632 G4int step(nEnergySteps);
633 while (step > 0)
634 {
635 step--;
636 G4double differentialCrossSection =
637 DifferentialCrossSection(particleDefinition,
638 k / eV,
639 value / eV,
640 shell);
641 if (differentialCrossSection >= crossSectionMaximum)
642 crossSectionMaximum = differentialCrossSection;
643 value *= stpEnergy;
644 }
645 //
646
647 G4double secondaryElectronKineticEnergy = 0.;
648 do
649 {
650 secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
651 } while(G4UniformRand()*crossSectionMaximum >
652 DifferentialCrossSection(particleDefinition, k/eV,
653 (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
654
655 return secondaryElectronKineticEnergy;
656
657 }
658
659 if (particleDefinition == G4Proton::ProtonDefinition())
660 {
661 G4double maximumKineticEnergyTransfer = 4.
662 * (electron_mass_c2 / proton_mass_c2) * k;
663
664 G4double crossSectionMaximum = 0.;
665 for (G4double value = waterStructure.IonisationEnergy(shell);
666 value <= 4. * waterStructure.IonisationEnergy(shell); value += 0.1 * eV)
667 {
668 G4double differentialCrossSection =
669 DifferentialCrossSection(particleDefinition,
670 k / eV,
671 value / eV,
672 shell);
673 if (differentialCrossSection >= crossSectionMaximum)
674 crossSectionMaximum = differentialCrossSection;
675 }
676
677 G4double secondaryElectronKineticEnergy = 0.;
678 do
679 {
680 secondaryElectronKineticEnergy = G4UniformRand()* maximumKineticEnergyTransfer;
681 } while(G4UniformRand()*crossSectionMaximum >=
682 DifferentialCrossSection(particleDefinition, k/eV,
683 (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
684
685 return secondaryElectronKineticEnergy;
686 }
687
688 return 0;
689}
690
691//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
693 G4double k,
694 G4double energyTransfer,
695 G4int ionizationLevelIndex)
696{
697 G4double sigma = 0.;
698
699 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)/eV)
700 {
701 G4double valueT1 = 0;
702 G4double valueT2 = 0;
703 G4double valueE21 = 0;
704 G4double valueE22 = 0;
705 G4double valueE12 = 0;
706 G4double valueE11 = 0;
707
708 G4double xs11 = 0;
709 G4double xs12 = 0;
710 G4double xs21 = 0;
711 G4double xs22 = 0;
712
713 if (particleDefinition == G4Electron::ElectronDefinition())
714 {
715
716 // Protection against out of boundary access - electron case : 1 MeV
717 if (k==eTdummyVec.back()) k=k*(1.-1e-12);
718 //
719
720 // k should be in eV and energy transfer eV also
721
722 auto t2 = std::upper_bound(eTdummyVec.begin(),
723 eTdummyVec.end(),
724 k);
725
726 auto t1 = t2 - 1;
727
728 // SI : the following condition avoids situations where energyTransfer >last vector element
729 if (energyTransfer <= eVecm[(*t1)].back()
730 && energyTransfer <= eVecm[(*t2)].back())
731 {
732 auto e12 =
733 std::upper_bound(eVecm[(*t1)].begin(),
734 eVecm[(*t1)].end(),
735 energyTransfer);
736 auto e11 = e12 - 1;
737
738 auto e22 =
739 std::upper_bound(eVecm[(*t2)].begin(),
740 eVecm[(*t2)].end(),
741 energyTransfer);
742 auto e21 = e22 - 1;
743
744 valueT1 = *t1;
745 valueT2 = *t2;
746 valueE21 = *e21;
747 valueE22 = *e22;
748 valueE12 = *e12;
749 valueE11 = *e11;
750
751 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
752 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
753 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
754 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
755
756 }
757
758 }
759
760 if (particleDefinition == G4Proton::ProtonDefinition())
761 {
762 // Protection against out of boundary access - proton case : 100 MeV
763 if (k==pTdummyVec.back()) k=k*(1.-1e-12);
764 //
765
766 // k should be in eV and energy transfer eV also
767 auto t2 = std::upper_bound(pTdummyVec.begin(),
768 pTdummyVec.end(),
769 k);
770 auto t1 = t2 - 1;
771
772 auto e12 = std::upper_bound(pVecm[(*t1)].begin(),
773 pVecm[(*t1)].end(),
774 energyTransfer);
775 auto e11 = e12 - 1;
776
777 auto e22 = std::upper_bound(pVecm[(*t2)].begin(),
778 pVecm[(*t2)].end(),
779 energyTransfer);
780 auto e21 = e22 - 1;
781
782 valueT1 = *t1;
783 valueT2 = *t2;
784 valueE21 = *e21;
785 valueE22 = *e22;
786 valueE12 = *e12;
787 valueE11 = *e11;
788
789 xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
790 xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
791 xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
792 xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
793
794 }
795
796 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
797 if (xsProduct != 0.)
798 {
799 sigma = QuadInterpolator(valueE11,
800 valueE12,
801 valueE21,
802 valueE22,
803 xs11,
804 xs12,
805 xs21,
806 xs22,
807 valueT1,
808 valueT2,
809 k,
810 energyTransfer);
811 }
812 }
813
814 return sigma;
815}
816
817//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
818
819G4double G4DNABornIonisationModel1::Interpolate(G4double e1,
820 G4double e2,
821 G4double e,
822 G4double xs1,
823 G4double xs2)
824{
825 G4double value = 0.;
826
827 // Log-log interpolation by default
828
829 if (e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
830 && !fasterCode)
831 {
832 G4double a = (std::log10(xs2) - std::log10(xs1))
833 / (std::log10(e2) - std::log10(e1));
834 G4double b = std::log10(xs2) - a * std::log10(e2);
835 G4double sigma = a * std::log10(e) + b;
836 value = (std::pow(10., sigma));
837 }
838
839 // Switch to lin-lin interpolation
840 /*
841 if ((e2-e1)!=0)
842 {
843 G4double d1 = xs1;
844 G4double d2 = xs2;
845 value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
846 }
847 */
848
849 // Switch to log-lin interpolation for faster code
850 if ((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
851 {
852 G4double d1 = std::log10(xs1);
853 G4double d2 = std::log10(xs2);
854 value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
855 }
856
857 // Switch to lin-lin interpolation for faster code
858 // in case one of xs1 or xs2 (=cum proba) value is zero
859
860 if ((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
861 {
862 G4double d1 = xs1;
863 G4double d2 = xs2;
864 value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
865 }
866
867 /*
868 G4cout
869 << e1 << " "
870 << e2 << " "
871 << e << " "
872 << xs1 << " "
873 << xs2 << " "
874 << value
875 << G4endl;
876 */
877
878 return value;
879}
880
881//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
882
883G4double G4DNABornIonisationModel1::QuadInterpolator(G4double e11,
884 G4double e12,
885 G4double e21,
886 G4double e22,
887 G4double xs11,
888 G4double xs12,
889 G4double xs21,
890 G4double xs22,
891 G4double t1,
892 G4double t2,
893 G4double t,
894 G4double e)
895{
896 G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
897 G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
898 G4double value = Interpolate(t1,
899 t2,
900 t,
901 interpolatedvalue1,
902 interpolatedvalue2);
903
904 return value;
905}
906
908 G4int level,
909 const G4ParticleDefinition* particle,
910 G4double kineticEnergy)
911{
912 std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
913 pos = tableData.find(particle->GetParticleName());
914
915 if (pos != tableData.end())
916 {
917 G4DNACrossSectionDataSet* table = pos->second;
918 return table->GetComponent(level)->FindValue(kineticEnergy);
919 }
920
921 return 0;
922}
923
924//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
925
926G4int G4DNABornIonisationModel1::RandomSelect(G4double k,
927 const G4String& particle)
928{
929 G4int level = 0;
930
931 std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
932 pos = tableData.find(particle);
933
934 if (pos != tableData.end())
935 {
936 G4DNACrossSectionDataSet* table = pos->second;
937
938 if (table != nullptr)
939 {
940 auto valuesBuffer = new G4double[table->NumberOfComponents()];
941 const auto n = (G4int)table->NumberOfComponents();
942 G4int i(n);
943 G4double value = 0.;
944
945 while (i > 0)
946 {
947 i--;
948 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
949 value += valuesBuffer[i];
950 }
951
952 value *= G4UniformRand();
953
954 i = n;
955
956 while (i > 0)
957 {
958 i--;
959
960 if (valuesBuffer[i] > value)
961 {
962 delete[] valuesBuffer;
963 return i;
964 }
965 value -= valuesBuffer[i];
966 }
967
968
969 delete[] valuesBuffer;
970
971 }
972 } else
973 {
974 G4Exception("G4DNABornIonisationModel1::RandomSelect",
975 "em0002",
977 "Model not applicable to particle type.");
978 }
979
980 return level;
981}
982
983//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
984
985G4double G4DNABornIonisationModel1::RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition* particleDefinition,
986 G4double k,
987 G4int shell)
988{
989 //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
990
991 G4double secondaryElectronKineticEnergy = 0.;
992
993 G4double random = G4UniformRand();
994
995 secondaryElectronKineticEnergy = TransferedEnergy(particleDefinition,
996 k / eV,
997 shell,
998 random) * eV
999 - waterStructure.IonisationEnergy(shell);
1000
1001 //G4cout << RandomTransferedEnergy(particleDefinition, k/eV, shell) << G4endl;
1002
1003 if (secondaryElectronKineticEnergy < 0.)
1004 return 0.;
1005 //
1006
1007 return secondaryElectronKineticEnergy;
1008}
1009
1010//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1011
1013 G4double k,
1014 G4int ionizationLevelIndex,
1015 G4double random)
1016{
1017 G4double nrj = 0.;
1018
1019 G4double valueK1 = 0;
1020 G4double valueK2 = 0;
1021 G4double valuePROB21 = 0;
1022 G4double valuePROB22 = 0;
1023 G4double valuePROB12 = 0;
1024 G4double valuePROB11 = 0;
1025
1026 G4double nrjTransf11 = 0;
1027 G4double nrjTransf12 = 0;
1028 G4double nrjTransf21 = 0;
1029 G4double nrjTransf22 = 0;
1030
1031 if (particleDefinition == G4Electron::ElectronDefinition())
1032 {
1033 // Protection against out of boundary access - electron case : 1 MeV
1034 if (k==eTdummyVec.back()) k=k*(1.-1e-12);
1035 //
1036
1037 // k should be in eV
1038 auto k2 = std::upper_bound(eTdummyVec.begin(),
1039 eTdummyVec.end(),
1040 k);
1041 auto k1 = k2 - 1;
1042
1043 /*
1044 G4cout << "----> k=" << k
1045 << " " << *k1
1046 << " " << *k2
1047 << " " << random
1048 << " " << ionizationLevelIndex
1049 << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
1050 << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
1051 << G4endl;
1052 */
1053
1054 // SI : the following condition avoids situations where random >last vector element
1055 if (random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
1056 && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back())
1057 {
1058 auto prob12 =
1059 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
1060 eProbaShellMap[ionizationLevelIndex][(*k1)].end(),
1061 random);
1062
1063 auto prob11 = prob12 - 1;
1064
1065 auto prob22 =
1066 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1067 eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1068 random);
1069
1070 auto prob21 = prob22 - 1;
1071
1072 valueK1 = *k1;
1073 valueK2 = *k2;
1074 valuePROB21 = *prob21;
1075 valuePROB22 = *prob22;
1076 valuePROB12 = *prob12;
1077 valuePROB11 = *prob11;
1078
1079 /*
1080 G4cout << " " << random << " " << valuePROB11 << " "
1081 << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1082 */
1083
1084 nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
1085 nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
1086 nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1087 nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1088
1089 /*
1090 G4cout << " " << ionizationLevelIndex << " "
1091 << random << " " <<valueK1 << " " << valueK2 << G4endl;
1092
1093 G4cout << " " << random << " " << nrjTransf11 << " "
1094 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1095 */
1096
1097 }
1098
1099 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1100 if (random > eProbaShellMap[ionizationLevelIndex][(*k1)].back())
1101 {
1102 auto prob22 =
1103 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1104 eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1105 random);
1106
1107 auto prob21 = prob22 - 1;
1108
1109 valueK1 = *k1;
1110 valueK2 = *k2;
1111 valuePROB21 = *prob21;
1112 valuePROB22 = *prob22;
1113
1114 //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1115
1116 nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1117 nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1118
1119 G4double interpolatedvalue2 = Interpolate(valuePROB21,
1120 valuePROB22,
1121 random,
1122 nrjTransf21,
1123 nrjTransf22);
1124
1125 // zeros are explicitly set
1126
1127 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1128
1129 /*
1130 G4cout << " " << ionizationLevelIndex << " "
1131 << random << " " <<valueK1 << " " << valueK2 << G4endl;
1132
1133 G4cout << " " << random << " " << nrjTransf11 << " "
1134 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1135
1136 G4cout << "ici" << " " << value << G4endl;
1137 */
1138
1139 return value;
1140 }
1141 }
1142 //
1143 else if (particleDefinition == G4Proton::ProtonDefinition())
1144 {
1145 // Protection against out of boundary access - proton case : 100 MeV
1146 if (k==pTdummyVec.back()) k=k*(1.-1e-12);
1147 //
1148
1149 // k should be in eV
1150
1151 auto k2 = std::upper_bound(pTdummyVec.begin(),
1152 pTdummyVec.end(),
1153 k);
1154
1155 auto k1 = k2 - 1;
1156
1157 /*
1158 G4cout << "----> k=" << k
1159 << " " << *k1
1160 << " " << *k2
1161 << " " << random
1162 << " " << ionizationLevelIndex
1163 << " " << pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1164 << " " << pProbaShellMap[ionizationLevelIndex][(*k2)].back()
1165 << G4endl;
1166 */
1167
1168 // SI : the following condition avoids situations where random > last vector element,
1169 // for eg. when the last element is zero
1170 if (random <= pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1171 && random <= pProbaShellMap[ionizationLevelIndex][(*k2)].back())
1172 {
1173 auto prob12 =
1174 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
1175 pProbaShellMap[ionizationLevelIndex][(*k1)].end(),
1176 random);
1177
1178 auto prob11 = prob12 - 1;
1179
1180 auto prob22 =
1181 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1182 pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1183 random);
1184
1185 auto prob21 = prob22 - 1;
1186
1187 valueK1 = *k1;
1188 valueK2 = *k2;
1189 valuePROB21 = *prob21;
1190 valuePROB22 = *prob22;
1191 valuePROB12 = *prob12;
1192 valuePROB11 = *prob11;
1193
1194 /*
1195 G4cout << " " << random << " " << valuePROB11 << " "
1196 << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1197 */
1198
1199 nrjTransf11 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
1200 nrjTransf12 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
1201 nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1202 nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1203
1204 /*
1205 G4cout << " " << ionizationLevelIndex << " "
1206 << random << " " <<valueK1 << " " << valueK2 << G4endl;
1207
1208 G4cout << " " << random << " " << nrjTransf11 << " "
1209 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1210 */
1211 }
1212
1213 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1214
1215 if (random > pProbaShellMap[ionizationLevelIndex][(*k1)].back())
1216 {
1217 auto prob22 =
1218 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1219 pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1220 random);
1221
1222 auto prob21 = prob22 - 1;
1223
1224 valueK1 = *k1;
1225 valueK2 = *k2;
1226 valuePROB21 = *prob21;
1227 valuePROB22 = *prob22;
1228
1229 //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1230
1231 nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1232 nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1233
1234 G4double interpolatedvalue2 = Interpolate(valuePROB21,
1235 valuePROB22,
1236 random,
1237 nrjTransf21,
1238 nrjTransf22);
1239
1240 // zeros are explicitly set
1241
1242 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1243
1244 /*
1245 G4cout << " " << ionizationLevelIndex << " "
1246 << random << " " <<valueK1 << " " << valueK2 << G4endl;
1247
1248 G4cout << " " << random << " " << nrjTransf11 << " "
1249 << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1250
1251 G4cout << "ici" << " " << value << G4endl;
1252 */
1253
1254 return value;
1255 }
1256 }
1257 // End electron and proton cases
1258
1259 G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21
1260 * nrjTransf22;
1261 //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
1262
1263 if (nrjTransfProduct != 0.)
1264 {
1265 nrj = QuadInterpolator(valuePROB11,
1266 valuePROB12,
1267 valuePROB21,
1268 valuePROB22,
1269 nrjTransf11,
1270 nrjTransf12,
1271 nrjTransf21,
1272 nrjTransf22,
1273 valueK1,
1274 valueK2,
1275 k,
1276 random);
1277 }
1278 //G4cout << nrj << endl;
1279
1280 return nrj;
1281}
@ eIonizedMolecule
G4double B(G4double temperature)
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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)
void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double) override
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void StartTracking(G4Track *) override
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNABornIonisationModel1(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")
static G4DNAChemistryManager * Instance()
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 DNAFast() const
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 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 *)