Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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);
71  fAtomDeexcitation = 0;
73  fpMolWaterDensity = 0;
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  {
215  SetLowEnergyLimit(lowEnergyLimit[electron]);
216  SetHighEnergyLimit(highEnergyLimit[electron]);
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
230  fpMolWaterDensity =
232  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
233 
234  //
235  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
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  //SI: additional protection if tcs interpolation method is modified
384  if (k<bindingEnergy) return;
385  //
386 
387  G4int Z = 8;
388  if(fAtomDeexcitation)
389  {
391 
392  if (ionizationShell <5 && ionizationShell >1)
393  {
394  as = G4AtomicShellEnumerator(4-ionizationShell);
395  }
396  else if (ionizationShell <2)
397  {
398  as = G4AtomicShellEnumerator(3);
399  }
400 
401  // FOR DEBUG ONLY
402  // if (ionizationShell == 4) {
403  //
404  // G4cout << "Z: " << Z << " as: " << as
405  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
406  // G4cout << "Press <Enter> key to continue..." << G4endl;
407  // G4cin.ignore();
408  // }
409 
410  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
411  secNumberInit = fvect->size();
412  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
413  secNumberFinal = fvect->size();
414  }
415 
416  G4double secondaryKinetic=-1000*eV;
417 
418  if (!fasterCode) secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
419 
420  // SI - 01/04/2014
421  if (fasterCode)
422  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
423  //
424 
425  G4ThreeVector deltaDirection =
426  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
427  Z, ionizationShell,
428  couple->GetMaterial());
429 
430  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
431 
432  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
433  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
434  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
435  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
436  finalPx /= finalMomentum;
437  finalPy /= finalMomentum;
438  finalPz /= finalMomentum;
439 
440  G4ThreeVector direction;
441  direction.set(finalPx,finalPy,finalPz);
442 
444 
445  // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
446  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
447  G4double deexSecEnergy = 0;
448  for (G4int j=secNumberInit; j < secNumberFinal; j++)
449  {
450  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
451  }
452 
453  if (!statCode)
454  {
456  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
457  }
458  else
459  {
462  }
463 
464  // SI - 01/04/2014
465  if (secondaryKinetic>0)
466  {
467  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
468  fvect->push_back(dp);
469  }
470  //
471 
472  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
474  ionizationShell,
475  theIncomingTrack);
476  }
477 
478 }
479 
480 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
481 
482 G4double
483 G4DNAEmfietzoglouIonisationModel::
484 RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
485  G4double k,
486  G4int shell)
487 {
488  // G4cout << "*** SLOW computation for "
489  // << " " << particleDefinition->GetParticleName() << G4endl;
490 
491  if(particleDefinition == G4Electron::ElectronDefinition())
492  {
493  G4double maximumEnergyTransfer = 0.;
494  if((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
495  maximumEnergyTransfer = k;
496  else
497  maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell))/ 2.;
498 
499  // SI : original method
500  /*
501  G4double crossSectionMaximum = 0.;
502  for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
503  {
504  G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
505  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
506  }
507  */
508 
509  // SI : alternative method
510  G4double crossSectionMaximum = 0.;
511 
512  G4double minEnergy = waterStructure.IonisationEnergy(shell);
513  G4double maxEnergy = maximumEnergyTransfer;
514  G4int nEnergySteps = 50;
515 
516  G4double value(minEnergy);
517  G4double stpEnergy(std::pow(maxEnergy / value,
518  1. / static_cast<G4double>(nEnergySteps - 1)));
519  G4int step(nEnergySteps);
520  while(step > 0)
521  {
522  step--;
523  G4double differentialCrossSection =
524  DifferentialCrossSection(particleDefinition,
525  k / eV,
526  value / eV,
527  shell);
528  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum =
529  differentialCrossSection;
530  value *= stpEnergy;
531  }
532  //
533 
534  G4double secondaryElectronKineticEnergy = 0.;
535  do
536  {
537  secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
538  }while(G4UniformRand()*crossSectionMaximum >
539  DifferentialCrossSection(particleDefinition, k/eV,
540  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
541 
542  return secondaryElectronKineticEnergy;
543 
544  }
545 
546  return 0;
547 }
548 
549 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
550 
551 // The following section is not used anymore but is kept for memory
552 // GetAngularDistribution()->SampleDirectionForShell is used instead
553 
554 /*
555  void G4DNAEmfietzoglouIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
556  G4double k,
557  G4double secKinetic,
558  G4double & cosTheta,
559  G4double & phi )
560  {
561  if (particleDefinition == G4Electron::ElectronDefinition())
562  {
563  phi = twopi * G4UniformRand();
564  if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
565  else if (secKinetic <= 200.*eV)
566  {
567  if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
568  else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
569  }
570  else
571  {
572  G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
573  cosTheta = std::sqrt(1.-sin2O);
574  }
575  }
576 
577  else if (particleDefinition == G4Proton::ProtonDefinition())
578  {
579  G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
580  phi = twopi * G4UniformRand();
581 
582  // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
583 
584  // Restriction below 100 eV from Emfietzoglou (2000)
585 
586  if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
587  else cosTheta = (2.*G4UniformRand())-1.;
588 
589  }
590  }
591  */
592 
593 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
595  G4double k,
596  G4double energyTransfer,
597  G4int ionizationLevelIndex)
598 {
599  G4double sigma = 0.;
600 
601  if(energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
602  {
603  G4double valueT1 = 0;
604  G4double valueT2 = 0;
605  G4double valueE21 = 0;
606  G4double valueE22 = 0;
607  G4double valueE12 = 0;
608  G4double valueE11 = 0;
609 
610  G4double xs11 = 0;
611  G4double xs12 = 0;
612  G4double xs21 = 0;
613  G4double xs22 = 0;
614 
615  if(particleDefinition == G4Electron::ElectronDefinition())
616  {
617  // k should be in eV and energy transfer eV also
618 
619  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
620  eTdummyVec.end(),
621  k);
622 
623  std::vector<double>::iterator t1 = t2 - 1;
624 
625  // SI : the following condition avoids situations where energyTransfer >last vector element
626  if(energyTransfer <= eVecm[(*t1)].back() && energyTransfer
627  <= eVecm[(*t2)].back())
628  {
629  std::vector<double>::iterator e12 =
630  std::upper_bound(eVecm[(*t1)].begin(),
631  eVecm[(*t1)].end(),
632  energyTransfer);
633  std::vector<double>::iterator e11 = e12 - 1;
634 
635  std::vector<double>::iterator e22 =
636  std::upper_bound(eVecm[(*t2)].begin(),
637  eVecm[(*t2)].end(),
638  energyTransfer);
639  std::vector<double>::iterator e21 = e22 - 1;
640 
641  valueT1 = *t1;
642  valueT2 = *t2;
643  valueE21 = *e21;
644  valueE22 = *e22;
645  valueE12 = *e12;
646  valueE11 = *e11;
647 
648  xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
649  xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
650  xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
651  xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
652 
653  //G4cout << "-------------------" << G4endl;
654  //G4cout << "ionizationLevelIndex=" << ionizationLevelIndex << G4endl;
655  //G4cout << "valueT1/eV=" << valueT1 << " valueT2/eV=" << valueT2 << G4endl;
656  //G4cout << "valueE11/eV=" << valueE11 << " valueE12/eV=" << valueE12 << " valueE21/eV=" << valueE21 << " valueE22/eV=" << valueE22 << G4endl;
657  //G4cout << "xs11=" << xs11 / ((1.e-22 / 3.343) * m*m) << G4endl;
658  //G4cout << "xs12=" << xs12 / ((1.e-22 / 3.343) * m*m) << G4endl;
659  //G4cout << "xs21=" << xs21 / ((1.e-22 / 3.343) * m*m) << G4endl;
660  //G4cout << "xs22=" << xs22 / ((1.e-22 / 3.343) * m*m) << G4endl;
661  //G4cout << "###################" << G4endl;
662 
663  }
664 
665  }
666 
667  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
668  if(xsProduct != 0.)
669  {
670  sigma = QuadInterpolator(valueE11,
671  valueE12,
672  valueE21,
673  valueE22,
674  xs11,
675  xs12,
676  xs21,
677  xs22,
678  valueT1,
679  valueT2,
680  k,
681  energyTransfer);
682  }
683 
684  }
685 
686  return sigma;
687 }
688 
689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
690 
691 G4double G4DNAEmfietzoglouIonisationModel::Interpolate(G4double e1,
692  G4double e2,
693  G4double e,
694  G4double xs1,
695  G4double xs2)
696 {
697 
698  G4double value = 0.;
699 
700  // Log-log interpolation by default
701 
702  if(e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
703  && !fasterCode)
704  {
705  G4double a = (std::log10(xs2) - std::log10(xs1))
706  / (std::log10(e2) - std::log10(e1));
707  G4double b = std::log10(xs2) - a * std::log10(e2);
708  G4double sigma = a * std::log10(e) + b;
709  value = (std::pow(10., sigma));
710  }
711 
712  // Switch to lin-lin interpolation
713  /*
714  if ((e2-e1)!=0)
715  {
716  G4double d1 = xs1;
717  G4double d2 = xs2;
718  value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
719  }
720  */
721 
722  // Switch to log-lin interpolation for faster code
723  if((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
724  {
725  G4double d1 = std::log10(xs1);
726  G4double d2 = std::log10(xs2);
727  value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
728  }
729 
730  // Switch to lin-lin interpolation for faster code
731  // in case one of xs1 or xs2 (=cum proba) value is zero
732 
733  if((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
734  {
735  G4double d1 = xs1;
736  G4double d2 = xs2;
737  value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
738  }
739 
740  /*
741  G4cout
742  << e1 << " "
743  << e2 << " "
744  << e << " "
745  << xs1 << " "
746  << xs2 << " "
747  << value
748  << G4endl;
749  */
750 
751  return value;
752 }
753 
754 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
755 
756 G4double G4DNAEmfietzoglouIonisationModel::QuadInterpolator(G4double e11,
757  G4double e12,
758  G4double e21,
759  G4double e22,
760  G4double xs11,
761  G4double xs12,
762  G4double xs21,
763  G4double xs22,
764  G4double t1,
765  G4double t2,
766  G4double t,
767  G4double e)
768 {
769  G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
770  G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
771  G4double value = Interpolate(t1,
772  t2,
773  t,
774  interpolatedvalue1,
775  interpolatedvalue2);
776 
777  return value;
778 }
779 
780 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
781 
782 G4int G4DNAEmfietzoglouIonisationModel::RandomSelect(G4double k,
783  const G4String& particle)
784 {
785  G4int level = 0;
786 
787  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
788  pos = tableData.find(particle);
789 
790  if(pos != tableData.end())
791  {
792  G4DNACrossSectionDataSet* table = pos->second;
793 
794  if(table != 0)
795  {
796  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
797  const size_t n(table->NumberOfComponents());
798  size_t i(n);
799  G4double value = 0.;
800 
801  while(i > 0)
802  {
803  i--;
804  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
805  value += valuesBuffer[i];
806  }
807 
808  value *= G4UniformRand();
809 
810  i = n;
811 
812  while(i > 0)
813  {
814  i--;
815 
816  if(valuesBuffer[i] > value)
817  {
818  delete[] valuesBuffer;
819  return i;
820  }
821  value -= valuesBuffer[i];
822  }
823 
824  if(valuesBuffer) delete[] valuesBuffer;
825 
826  }
827  }
828  else
829  {
830  G4Exception("G4DNAEmfietzoglouIonisationModel::RandomSelect",
831  "em0002",
833  "Model not applicable to particle type.");
834  }
835 
836  return level;
837 }
838 
839 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
840 
841 G4double G4DNAEmfietzoglouIonisationModel::RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition* particleDefinition,
842  G4double k,
843  G4int shell)
844 {
845  //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
846 
847  G4double secondaryElectronKineticEnergy = 0.;
848 
849  secondaryElectronKineticEnergy = RandomTransferedEnergy(particleDefinition,
850  k / eV,
851  shell)
852  * eV
853  - waterStructure.IonisationEnergy(shell);
854 
855  //G4cout << RandomTransferedEnergy(particleDefinition, k/eV, shell) << G4endl;
856  // SI - 01/04/2014
857  if(secondaryElectronKineticEnergy < 0.) return 0.;
858  //
859 
860  return secondaryElectronKineticEnergy;
861 }
862 
863 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
864 
865 G4double G4DNAEmfietzoglouIonisationModel::RandomTransferedEnergy(G4ParticleDefinition* particleDefinition,
866  G4double k,
867  G4int ionizationLevelIndex)
868 {
869 
870  G4double random = G4UniformRand();
871 
872  G4double nrj = 0.;
873 
874  G4double valueK1 = 0;
875  G4double valueK2 = 0;
876  G4double valuePROB21 = 0;
877  G4double valuePROB22 = 0;
878  G4double valuePROB12 = 0;
879  G4double valuePROB11 = 0;
880 
881  G4double nrjTransf11 = 0;
882  G4double nrjTransf12 = 0;
883  G4double nrjTransf21 = 0;
884  G4double nrjTransf22 = 0;
885 
886  if (particleDefinition == G4Electron::ElectronDefinition())
887  {
888  // k should be in eV
889  std::vector<double>::iterator k2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
890 
891  std::vector<double>::iterator k1 = k2-1;
892 
893  /*
894  G4cout << "----> k=" << k
895  << " " << *k1
896  << " " << *k2
897  << " " << random
898  << " " << ionizationLevelIndex
899  << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
900  << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
901  << G4endl;
902  */
903 
904  // SI : the following condition avoids situations where random >last vector element
905  if ( random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
906  && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back() )
907 
908  {
909  std::vector<double>::iterator prob12 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
910  eProbaShellMap[ionizationLevelIndex][(*k1)].end(), random);
911 
912  std::vector<double>::iterator prob11 = prob12-1;
913 
914  std::vector<double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
915  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
916 
917  std::vector<double>::iterator prob21 = prob22-1;
918 
919  valueK1 =*k1;
920  valueK2 =*k2;
921  valuePROB21 =*prob21;
922  valuePROB22 =*prob22;
923  valuePROB12 =*prob12;
924  valuePROB11 =*prob11;
925 
926  /*
927  G4cout << " " << random << " " << valuePROB11 << " "
928  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
929  */
930 
931  nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
932  nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
933  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
934  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
935 
936  /*
937  G4cout << " " << ionizationLevelIndex << " "
938  << random << " " <<valueK1 << " " << valueK2 << G4endl;
939 
940  G4cout << " " << random << " " << nrjTransf11 << " "
941  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
942  */
943 
944  }
945 
946  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
947 
948  if ( random > eProbaShellMap[ionizationLevelIndex][(*k1)].back() )
949 
950  {
951  std::vector<double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
952  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
953 
954  std::vector<double>::iterator prob21 = prob22-1;
955 
956  valueK1 =*k1;
957  valueK2 =*k2;
958  valuePROB21 =*prob21;
959  valuePROB22 =*prob22;
960 
961  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
962 
963  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
964  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
965 
966  G4double interpolatedvalue2 = Interpolate(valuePROB21, valuePROB22, random, nrjTransf21, nrjTransf22);
967 
968  // zeros are explicitely set
969 
970  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
971 
972  /*
973  G4cout << " " << ionizationLevelIndex << " "
974  << random << " " <<valueK1 << " " << valueK2 << G4endl;
975 
976  G4cout << " " << random << " " << nrjTransf11 << " "
977  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
978 
979  G4cout << "ici" << " " << value << G4endl;
980  */
981 
982  return value;
983  }
984 
985  }
986 
987  // End electron
988 
989  G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21 * nrjTransf22;
990 
991  //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
992 
993  if (nrjTransfProduct != 0.)
994  {
995  nrj = QuadInterpolator( valuePROB11, valuePROB12,
996  valuePROB21, valuePROB22,
997  nrjTransf11, nrjTransf12,
998  nrjTransf21, nrjTransf22,
999  valueK1, valueK2,
1000  k, random);
1001  }
1002 
1003  //G4cout << nrj << endl;
1004 
1005  return nrj;
1006 }
void set(double x, double y, double z)
G4DNAEmfietzoglouIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static G4LossTableManager * Instance()
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
double x() const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
size_t GetIndex() const
Definition: G4Material.hh:262
static const G4double d2
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
virtual G4bool LoadData(const G4String &argFileName)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double electron_mass_c2
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
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
#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)
static constexpr double m
Definition: G4SIunits.hh:129
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
static const G4double d1
virtual size_t NumberOfComponents(void) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static G4DNAMolecularMaterial * Instance()
G4double GetPDGMass() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
Hep3Vector unit() const
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
double y() const
const G4Track * GetCurrentTrack() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
G4double bindingEnergy(G4int A, G4int Z)
static constexpr double keV
Definition: G4SIunits.hh:216
static const G4double pos
G4AtomicShellEnumerator
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const