Geant4  10.02
G4DNAEmfietzoglouIonisationModel.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 // Based on the work described in
27 // Rad Res 163, 98-111 (2005)
28 // D. Emfietzoglou, H. Nikjoo
29 //
30 // Authors of the class (2014):
31 // I. Kyriakou (kyriak@cc.uoi.gr)
32 // D. Emfietzoglou (demfietz@cc.uoi.gr)
33 // S. Incerti (incerti@cenbg.in2p3.fr)
34 //
35 
37 #include "G4PhysicalConstants.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4UAtomicDeexcitation.hh"
40 #include "G4LossTableManager.hh"
41 #include "G4DNAChemistryManager.hh"
43 #include "G4DNABornAngle.hh"
44 #include "G4DeltaAngle.hh"
45 
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 using namespace std;
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51 
53  const G4String& nam) :
54  G4VEmModel(nam), isInitialised(false)
55 {
56  verboseLevel = 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, sampling of atoms
62  // 4 = entering in methods
63 
64  if(verboseLevel > 0)
65  {
66  G4cout << "Emfietzoglou ionisation model is constructed " << G4endl;
67  }
68 
69  //Mark this model as "applicable" for atomic deexcitation
70  SetDeexcitationFlag(true);
74 
75  // define default angular generator
77 
78  SetLowEnergyLimit(10. * eV);
79  SetHighEnergyLimit(10. * keV);
80 
81  // Selection of computation method
82 
83  fasterCode = false;
84 
85  // Selection of stationary mode
86 
87  statCode = false;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
94  // Cross section
95 
96  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
97  for(pos = tableData.begin(); pos != tableData.end(); ++pos)
98  {
99  G4DNACrossSectionDataSet* table = pos->second;
100  delete table;
101  }
102 
103  // Final state
104 
105  eVecm.clear();
106 
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110 
112  const G4DataVector& /*cuts*/)
113 {
114 
115  if(verboseLevel > 3)
116  {
117  G4cout << "Calling G4DNAEmfietzoglouIonisationModel::Initialise()" << G4endl;
118  }
119 
120  // Energy limits
121 
122  G4String fileElectron("dna/sigma_ionisation_e_emfietzoglou");
123 
125 
127 
128  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
129 
130  char *path = getenv("G4LEDATA");
131 
132  // *** ELECTRON
133 
134  electron = electronDef->GetParticleName();
135 
136  tableFile[electron] = fileElectron;
137 
138  lowEnergyLimit[electron] = 10. * eV;
139  highEnergyLimit[electron] = 10. * keV;
140 
141  // Cross section
142 
144  tableE->LoadData(fileElectron);
145 
146  tableData[electron] = tableE;
147 
148  // Final state
149 
150  std::ostringstream eFullFileName;
151 
152  if (fasterCode) eFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_e_emfietzoglou.dat";
153  if (!fasterCode) eFullFileName << path << "/dna/sigmadiff_ionisation_e_emfietzoglou.dat";
154 
155  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
156 
157  if (!eDiffCrossSection)
158  {
159  if (fasterCode) G4Exception("G4DNAEmfietzoglouIonisationModel::Initialise","em0003",
160  FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_e_emfietzoglou.dat");
161 
162  if (!fasterCode) G4Exception("G4DNAEmfietzoglouIonisationModel::Initialise","em0003",
163  FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_emfietzoglou.dat");
164  }
165 
166  //
167 
168  // Clear the arrays for re-initialization case (MT mode)
169  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
170 
171  eTdummyVec.clear();
172 
173  eVecm.clear();
174 
175  eProbaShellMap->clear();
176 
177  eDiffCrossSectionData->clear();
178 
179  eNrjTransfData->clear();
180 
181  //
182 
183  eTdummyVec.push_back(0.);
184  while(!eDiffCrossSection.eof())
185  {
186  double tDummy;
187  double eDummy;
188  eDiffCrossSection>>tDummy>>eDummy;
189  if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
190  for (int j=0; j<5; j++)
191  {
192  eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
193 
194  if (fasterCode)
195  {
196  eNrjTransfData[j][tDummy][eDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
197  eProbaShellMap[j][tDummy].push_back(eDiffCrossSectionData[j][tDummy][eDummy]);
198  }
199 
200  // SI - only if eof is not reached
201  if (!eDiffCrossSection.eof() && !fasterCode)
202  {
203  eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
204  }
205 
206  if (!fasterCode) eVecm[tDummy].push_back(eDummy);
207 
208  }
209  }
210 
211  //
212 
213  if (particle==electronDef)
214  {
217  }
218 
219  if( verboseLevel>0 )
220  {
221  G4cout << "Emfietzoglou ionisation model is initialized " << G4endl
222  << "Energy range: "
223  << LowEnergyLimit() / eV << " eV - "
224  << HighEnergyLimit() / keV << " keV for "
225  << particle->GetParticleName()
226  << G4endl;
227  }
228 
229  // Initialize water density pointer
232  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
233 
234  //
236 
237  if (isInitialised)
238  { return;}
240  isInitialised = true;
241 }
242 
243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244 
247  const G4ParticleDefinition* particleDefinition,
248  G4double ekin,
249  G4double,
250  G4double)
251 {
252  if(verboseLevel > 3)
253  {
254  G4cout
255  << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouIonisationModel"
256  << G4endl;
257  }
258 
259  if (particleDefinition != G4Electron::ElectronDefinition()) return 0; // necessary ??
260 
261  // Calculate total cross section for model
262 
263  G4double lowLim = 0;
264  G4double highLim = 0;
265  G4double sigma=0;
266 
267  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
268 
269  if(waterDensity!= 0.0)
270  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
271  {
272  const G4String& particleName = particleDefinition->GetParticleName();
273 
274  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
275  pos1 = lowEnergyLimit.find(particleName);
276  if (pos1 != lowEnergyLimit.end())
277  {
278  lowLim = pos1->second;
279  }
280 
281  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
282  pos2 = highEnergyLimit.find(particleName);
283  if (pos2 != highEnergyLimit.end())
284  {
285  highLim = pos2->second;
286  }
287 
288  if (ekin >= lowLim && ekin < highLim)
289  {
290  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
291  pos = tableData.find(particleName);
292 
293  if (pos != tableData.end())
294  {
295  G4DNACrossSectionDataSet* table = pos->second;
296  if (table != 0)
297  {
298  sigma = table->FindValue(ekin);
299  }
300  }
301  else
302  {
303  G4Exception("G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume","em0002",
304  FatalException,"Model not applicable to particle type.");
305  }
306  }
307 
308  if (verboseLevel > 2)
309  {
310  G4cout << "__________________________________" << G4endl;
311  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO START" << G4endl;
312  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
313  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
314  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
315  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO END" << G4endl;
316  }
317  } // if (waterMaterial)
318 
319  return sigma*waterDensity;
320 }
321 
322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
323 
325 SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
326  const G4MaterialCutsCouple* couple,
327  const G4DynamicParticle* particle,
328  G4double,
329  G4double)
330 {
331 
332  if(verboseLevel > 3)
333  {
334  G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouIonisationModel"
335  << G4endl;
336  }
337 
338  G4double lowLim = 0;
339  G4double highLim = 0;
340 
341  G4double k = particle->GetKineticEnergy();
342 
343  const G4String& particleName = particle->GetDefinition()->GetParticleName();
344 
345  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
346  pos1 = lowEnergyLimit.find(particleName);
347 
348  if (pos1 != lowEnergyLimit.end())
349  {
350  lowLim = pos1->second;
351  }
352 
353  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
354  pos2 = highEnergyLimit.find(particleName);
355 
356  if (pos2 != highEnergyLimit.end())
357  {
358  highLim = pos2->second;
359  }
360 
361  if (k >= lowLim && k < highLim)
362  {
363  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
364  G4double particleMass = particle->GetDefinition()->GetPDGMass();
365  G4double totalEnergy = k + particleMass;
366  G4double pSquare = k * (totalEnergy + particleMass);
367  G4double totalMomentum = std::sqrt(pSquare);
368 
369  G4int ionizationShell = 0;
370 
371  ionizationShell = RandomSelect(k,particleName);
372 
373  // AM: sample deexcitation
374  // here we assume that H_{2}O electronic levels are the same as Oxygen.
375  // this can be considered true with a rough 10% error in energy on K-shell,
376 
377  G4int secNumberInit = 0;// need to know at a certain point the energy of secondaries
378  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
379 
381  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
382 
383  G4int Z = 8;
385  {
387 
388  if (ionizationShell <5 && ionizationShell >1)
389  {
390  as = G4AtomicShellEnumerator(4-ionizationShell);
391  }
392  else if (ionizationShell <2)
393  {
394  as = G4AtomicShellEnumerator(3);
395  }
396 
397  // FOR DEBUG ONLY
398  // if (ionizationShell == 4) {
399  //
400  // G4cout << "Z: " << Z << " as: " << as
401  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
402  // G4cout << "Press <Enter> key to continue..." << G4endl;
403  // G4cin.ignore();
404  // }
405 
406  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
407  secNumberInit = fvect->size();
408  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
409  secNumberFinal = fvect->size();
410  }
411 
412  G4double secondaryKinetic=-1000*eV;
413 
414  if (!fasterCode) secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
415 
416  // SI - 01/04/2014
417  if (fasterCode)
418  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
419  //
420 
421  G4ThreeVector deltaDirection =
422  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
423  Z, ionizationShell,
424  couple->GetMaterial());
425 
426  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
427 
428  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
429  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
430  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
431  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
432  finalPx /= finalMomentum;
433  finalPy /= finalMomentum;
434  finalPz /= finalMomentum;
435 
436  G4ThreeVector direction;
437  direction.set(finalPx,finalPy,finalPz);
438 
440 
441  // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
442  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
443  G4double deexSecEnergy = 0;
444  for (G4int j=secNumberInit; j < secNumberFinal; j++)
445  {
446  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
447  }
448 
449  if (!statCode)
450  {
452  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
453  }
454  else
455  {
458  }
459 
460  // SI - 01/04/2014
461  if (secondaryKinetic>0)
462  {
463  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
464  fvect->push_back(dp);
465  }
466  //
467 
468  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
470  ionizationShell,
471  theIncomingTrack);
472  }
473 
474 }
475 
476 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
477 
478 G4double
481  G4double k,
482  G4int shell)
483 {
484  // G4cout << "*** SLOW computation for "
485  // << " " << particleDefinition->GetParticleName() << G4endl;
486 
487  if(particleDefinition == G4Electron::ElectronDefinition())
488  {
489  G4double maximumEnergyTransfer = 0.;
490  if((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
491  maximumEnergyTransfer = k;
492  else
493  maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell))/ 2.;
494 
495  // SI : original method
496  /*
497  G4double crossSectionMaximum = 0.;
498  for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
499  {
500  G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
501  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
502  }
503  */
504 
505  // SI : alternative method
506  G4double crossSectionMaximum = 0.;
507 
508  G4double minEnergy = waterStructure.IonisationEnergy(shell);
509  G4double maxEnergy = maximumEnergyTransfer;
510  G4int nEnergySteps = 50;
511 
512  G4double value(minEnergy);
513  G4double stpEnergy(std::pow(maxEnergy / value,
514  1. / static_cast<G4double>(nEnergySteps - 1)));
515  G4int step(nEnergySteps);
516  while(step > 0)
517  {
518  step--;
519  G4double differentialCrossSection =
520  DifferentialCrossSection(particleDefinition,
521  k / eV,
522  value / eV,
523  shell);
524  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum =
525  differentialCrossSection;
526  value *= stpEnergy;
527  }
528  //
529 
530  G4double secondaryElectronKineticEnergy = 0.;
531  do
532  {
533  secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
534  }while(G4UniformRand()*crossSectionMaximum >
535  DifferentialCrossSection(particleDefinition, k/eV,
536  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
537 
538  return secondaryElectronKineticEnergy;
539 
540  }
541 
542  return 0;
543 }
544 
545 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
546 
547 // The following section is not used anymore but is kept for memory
548 // GetAngularDistribution()->SampleDirectionForShell is used instead
549 
550 /*
551  void G4DNAEmfietzoglouIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
552  G4double k,
553  G4double secKinetic,
554  G4double & cosTheta,
555  G4double & phi )
556  {
557  if (particleDefinition == G4Electron::ElectronDefinition())
558  {
559  phi = twopi * G4UniformRand();
560  if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
561  else if (secKinetic <= 200.*eV)
562  {
563  if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
564  else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
565  }
566  else
567  {
568  G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
569  cosTheta = std::sqrt(1.-sin2O);
570  }
571  }
572 
573  else if (particleDefinition == G4Proton::ProtonDefinition())
574  {
575  G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
576  phi = twopi * G4UniformRand();
577 
578  // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
579 
580  // Restriction below 100 eV from Emfietzoglou (2000)
581 
582  if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
583  else cosTheta = (2.*G4UniformRand())-1.;
584 
585  }
586  }
587  */
588 
589 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
591  G4double k,
592  G4double energyTransfer,
593  G4int ionizationLevelIndex)
594 {
595  G4double sigma = 0.;
596 
597  if(energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
598  {
599  G4double valueT1 = 0;
600  G4double valueT2 = 0;
601  G4double valueE21 = 0;
602  G4double valueE22 = 0;
603  G4double valueE12 = 0;
604  G4double valueE11 = 0;
605 
606  G4double xs11 = 0;
607  G4double xs12 = 0;
608  G4double xs21 = 0;
609  G4double xs22 = 0;
610 
611  if(particleDefinition == G4Electron::ElectronDefinition())
612  {
613  // k should be in eV and energy transfer eV also
614 
615  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
616  eTdummyVec.end(),
617  k);
618 
619  std::vector<double>::iterator t1 = t2 - 1;
620 
621  // SI : the following condition avoids situations where energyTransfer >last vector element
622  if(energyTransfer <= eVecm[(*t1)].back() && energyTransfer
623  <= eVecm[(*t2)].back())
624  {
625  std::vector<double>::iterator e12 =
626  std::upper_bound(eVecm[(*t1)].begin(),
627  eVecm[(*t1)].end(),
628  energyTransfer);
629  std::vector<double>::iterator e11 = e12 - 1;
630 
631  std::vector<double>::iterator e22 =
632  std::upper_bound(eVecm[(*t2)].begin(),
633  eVecm[(*t2)].end(),
634  energyTransfer);
635  std::vector<double>::iterator e21 = e22 - 1;
636 
637  valueT1 = *t1;
638  valueT2 = *t2;
639  valueE21 = *e21;
640  valueE22 = *e22;
641  valueE12 = *e12;
642  valueE11 = *e11;
643 
644  xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
645  xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
646  xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
647  xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
648 
649  //G4cout << "-------------------" << G4endl;
650  //G4cout << "ionizationLevelIndex=" << ionizationLevelIndex << G4endl;
651  //G4cout << "valueT1/eV=" << valueT1 << " valueT2/eV=" << valueT2 << G4endl;
652  //G4cout << "valueE11/eV=" << valueE11 << " valueE12/eV=" << valueE12 << " valueE21/eV=" << valueE21 << " valueE22/eV=" << valueE22 << G4endl;
653  //G4cout << "xs11=" << xs11 / ((1.e-22 / 3.343) * m*m) << G4endl;
654  //G4cout << "xs12=" << xs12 / ((1.e-22 / 3.343) * m*m) << G4endl;
655  //G4cout << "xs21=" << xs21 / ((1.e-22 / 3.343) * m*m) << G4endl;
656  //G4cout << "xs22=" << xs22 / ((1.e-22 / 3.343) * m*m) << G4endl;
657  //G4cout << "###################" << G4endl;
658 
659  }
660 
661  }
662 
663  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
664  if(xsProduct != 0.)
665  {
666  sigma = QuadInterpolator(valueE11,
667  valueE12,
668  valueE21,
669  valueE22,
670  xs11,
671  xs12,
672  xs21,
673  xs22,
674  valueT1,
675  valueT2,
676  k,
677  energyTransfer);
678  }
679 
680  }
681 
682  return sigma;
683 }
684 
685 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
686 
688  G4double e2,
689  G4double e,
690  G4double xs1,
691  G4double xs2)
692 {
693 
694  G4double value = 0.;
695 
696  // Log-log interpolation by default
697 
698  if(e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
699  && !fasterCode)
700  {
701  G4double a = (std::log10(xs2) - std::log10(xs1))
702  / (std::log10(e2) - std::log10(e1));
703  G4double b = std::log10(xs2) - a * std::log10(e2);
704  G4double sigma = a * std::log10(e) + b;
705  value = (std::pow(10., sigma));
706  }
707 
708  // Switch to lin-lin interpolation
709  /*
710  if ((e2-e1)!=0)
711  {
712  G4double d1 = xs1;
713  G4double d2 = xs2;
714  value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
715  }
716  */
717 
718  // Switch to log-lin interpolation for faster code
719  if((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
720  {
721  G4double d1 = std::log10(xs1);
722  G4double d2 = std::log10(xs2);
723  value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
724  }
725 
726  // Switch to lin-lin interpolation for faster code
727  // in case one of xs1 or xs2 (=cum proba) value is zero
728 
729  if((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
730  {
731  G4double d1 = xs1;
732  G4double d2 = xs2;
733  value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
734  }
735 
736  /*
737  G4cout
738  << e1 << " "
739  << e2 << " "
740  << e << " "
741  << xs1 << " "
742  << xs2 << " "
743  << value
744  << G4endl;
745  */
746 
747  return value;
748 }
749 
750 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
751 
753  G4double e12,
754  G4double e21,
755  G4double e22,
756  G4double xs11,
757  G4double xs12,
758  G4double xs21,
759  G4double xs22,
760  G4double t1,
761  G4double t2,
762  G4double t,
763  G4double e)
764 {
765  G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
766  G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
767  G4double value = Interpolate(t1,
768  t2,
769  t,
770  interpolatedvalue1,
771  interpolatedvalue2);
772 
773  return value;
774 }
775 
776 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
777 
779  const G4String& particle)
780 {
781  G4int level = 0;
782 
783  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
784  pos = tableData.find(particle);
785 
786  if(pos != tableData.end())
787  {
788  G4DNACrossSectionDataSet* table = pos->second;
789 
790  if(table != 0)
791  {
792  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
793  const size_t n(table->NumberOfComponents());
794  size_t i(n);
795  G4double value = 0.;
796 
797  while(i > 0)
798  {
799  i--;
800  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
801  value += valuesBuffer[i];
802  }
803 
804  value *= G4UniformRand();
805 
806  i = n;
807 
808  while(i > 0)
809  {
810  i--;
811 
812  if(valuesBuffer[i] > value)
813  {
814  delete[] valuesBuffer;
815  return i;
816  }
817  value -= valuesBuffer[i];
818  }
819 
820  if(valuesBuffer) delete[] valuesBuffer;
821 
822  }
823  }
824  else
825  {
826  G4Exception("G4DNAEmfietzoglouIonisationModel::RandomSelect",
827  "em0002",
829  "Model not applicable to particle type.");
830  }
831 
832  return level;
833 }
834 
835 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
836 
838  G4double k,
839  G4int shell)
840 {
841  //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
842 
843  G4double secondaryElectronKineticEnergy = 0.;
844 
845  secondaryElectronKineticEnergy = RandomTransferedEnergy(particleDefinition,
846  k / eV,
847  shell)
848  * eV
850 
851  //G4cout << RandomTransferedEnergy(particleDefinition, k/eV, shell) << G4endl;
852  // SI - 01/04/2014
853  if(secondaryElectronKineticEnergy < 0.) return 0.;
854  //
855 
856  return secondaryElectronKineticEnergy;
857 }
858 
859 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
860 
862  G4double k,
863  G4int ionizationLevelIndex)
864 {
865 
866  G4double random = G4UniformRand();
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  if (particleDefinition == G4Electron::ElectronDefinition())
883  {
884  // k should be in eV
885  std::vector<double>::iterator k2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
886 
887  std::vector<double>::iterator k1 = k2-1;
888 
889  /*
890  G4cout << "----> k=" << k
891  << " " << *k1
892  << " " << *k2
893  << " " << random
894  << " " << ionizationLevelIndex
895  << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
896  << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
897  << G4endl;
898  */
899 
900  // SI : the following condition avoids situations where random >last vector element
901  if ( random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
902  && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back() )
903 
904  {
905  std::vector<double>::iterator prob12 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
906  eProbaShellMap[ionizationLevelIndex][(*k1)].end(), random);
907 
908  std::vector<double>::iterator prob11 = prob12-1;
909 
910  std::vector<double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
911  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
912 
913  std::vector<double>::iterator prob21 = prob22-1;
914 
915  valueK1 =*k1;
916  valueK2 =*k2;
917  valuePROB21 =*prob21;
918  valuePROB22 =*prob22;
919  valuePROB12 =*prob12;
920  valuePROB11 =*prob11;
921 
922  /*
923  G4cout << " " << random << " " << valuePROB11 << " "
924  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
925  */
926 
927  nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
928  nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
929  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
930  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
931 
932  /*
933  G4cout << " " << ionizationLevelIndex << " "
934  << random << " " <<valueK1 << " " << valueK2 << G4endl;
935 
936  G4cout << " " << random << " " << nrjTransf11 << " "
937  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
938  */
939 
940  }
941 
942  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
943 
944  if ( random > eProbaShellMap[ionizationLevelIndex][(*k1)].back() )
945 
946  {
947  std::vector<double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
948  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
949 
950  std::vector<double>::iterator prob21 = prob22-1;
951 
952  valueK1 =*k1;
953  valueK2 =*k2;
954  valuePROB21 =*prob21;
955  valuePROB22 =*prob22;
956 
957  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
958 
959  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
960  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
961 
962  G4double interpolatedvalue2 = Interpolate(valuePROB21, valuePROB22, random, nrjTransf21, nrjTransf22);
963 
964  // zeros are explicitely set
965 
966  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
967 
968  /*
969  G4cout << " " << ionizationLevelIndex << " "
970  << random << " " <<valueK1 << " " << valueK2 << G4endl;
971 
972  G4cout << " " << random << " " << nrjTransf11 << " "
973  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
974 
975  G4cout << "ici" << " " << value << G4endl;
976  */
977 
978  return value;
979  }
980 
981  }
982 
983  // End electron
984 
985  G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21 * nrjTransf22;
986 
987  //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
988 
989  if (nrjTransfProduct != 0.)
990  {
991  nrj = QuadInterpolator( valuePROB11, valuePROB12,
992  valuePROB21, valuePROB22,
993  nrjTransf11, nrjTransf12,
994  nrjTransf21, nrjTransf22,
995  valueK1, valueK2,
996  k, random);
997  }
998 
999  //G4cout << nrj << endl;
1000 
1001  return nrj;
1002 }
static const double cm
Definition: G4SIunits.hh:118
G4DNAEmfietzoglouIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4DNAEmfietzoglouWaterIonisationStructure waterStructure
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
static const G4double d1
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static G4LossTableManager * Instance()
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
static const G4double e2
virtual G4bool LoadData(const G4String &argFileName)
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4ParticleDefinition * GetDefinition() const
G4double a
Definition: TRTMaterials.hh:39
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const std::vector< G4double > * fpMolWaterDensity
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
const G4int n
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4double RandomTransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
G4int RandomSelect(G4double energy, const G4String &particle)
virtual size_t NumberOfComponents(void) const
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
G4double GetPDGMass() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
Method used by DNA physics model to create a water molecule.
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
const G4Track * GetCurrentTrack() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
G4VAtomDeexcitation * AtomDeexcitation()
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
static const G4double d2
G4ThreeVector G4ParticleMomentum
G4double bindingEnergy(G4int A, G4int Z)
static const G4double pos
G4AtomicShellEnumerator
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
const G4Material * GetMaterial() const