Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNAEmfietzoglouIonisationModel.cc
이 파일의 문서화 페이지로 가기
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  G4double tDummy;
187  G4double eDummy;
188  eDiffCrossSection>>tDummy>>eDummy;
189  if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
190  for (G4int 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
230 
233  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
234 
235  // AD
236 
238 
239  if (isInitialised)
240  { return;}
242  isInitialised = true;
243 }
244 
245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246 
249  const G4ParticleDefinition* particleDefinition,
250  G4double ekin,
251  G4double,
252  G4double)
253 {
254  if(verboseLevel > 3)
255  {
256  G4cout
257  << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouIonisationModel"
258  << G4endl;
259  }
260 
261  if (particleDefinition != G4Electron::ElectronDefinition()) return 0; // necessary ??
262 
263  // Calculate total cross section for model
264 
265  G4double lowLim = 0;
266  G4double highLim = 0;
267  G4double sigma=0;
268 
269  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
270 
271  if(waterDensity!= 0.0)
272  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
273  {
274  const G4String& particleName = particleDefinition->GetParticleName();
275 
276  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
277  pos1 = lowEnergyLimit.find(particleName);
278  if (pos1 != lowEnergyLimit.end())
279  {
280  lowLim = pos1->second;
281  }
282 
283  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
284  pos2 = highEnergyLimit.find(particleName);
285  if (pos2 != highEnergyLimit.end())
286  {
287  highLim = pos2->second;
288  }
289 
290  if (ekin >= lowLim && ekin < highLim)
291  {
292  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
293  pos = tableData.find(particleName);
294 
295  if (pos != tableData.end())
296  {
297  G4DNACrossSectionDataSet* table = pos->second;
298  if (table != 0)
299  {
300  sigma = table->FindValue(ekin);
301  }
302  }
303  else
304  {
305  G4Exception("G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume","em0002",
306  FatalException,"Model not applicable to particle type.");
307  }
308  }
309 
310  if (verboseLevel > 2)
311  {
312  G4cout << "__________________________________" << G4endl;
313  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO START" << G4endl;
314  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
315  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
316  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
317  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO END" << G4endl;
318  }
319  } // if (waterMaterial)
320 
321  return sigma*waterDensity;
322 }
323 
324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325 
327 SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
328  const G4MaterialCutsCouple* couple,
329  const G4DynamicParticle* particle,
330  G4double,
331  G4double)
332 {
333 
334  if(verboseLevel > 3)
335  {
336  G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouIonisationModel"
337  << G4endl;
338  }
339 
340  G4double lowLim = 0;
341  G4double highLim = 0;
342 
343  G4double k = particle->GetKineticEnergy();
344 
345  const G4String& particleName = particle->GetDefinition()->GetParticleName();
346 
347  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
348  pos1 = lowEnergyLimit.find(particleName);
349 
350  if (pos1 != lowEnergyLimit.end())
351  {
352  lowLim = pos1->second;
353  }
354 
355  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
356  pos2 = highEnergyLimit.find(particleName);
357 
358  if (pos2 != highEnergyLimit.end())
359  {
360  highLim = pos2->second;
361  }
362 
363  if (k >= lowLim && k < highLim)
364  {
365  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
366  G4double particleMass = particle->GetDefinition()->GetPDGMass();
367  G4double totalEnergy = k + particleMass;
368  G4double pSquare = k * (totalEnergy + particleMass);
369  G4double totalMomentum = std::sqrt(pSquare);
370 
371  G4int ionizationShell = 0;
372 
373  ionizationShell = RandomSelect(k,particleName);
374 
376  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
377 
378  // SI : additional protection if tcs interpolation method is modified
379  if (k<bindingEnergy) return;
380  //
381 
382  G4double secondaryKinetic=-1000*eV;
383 
384  if (!fasterCode) secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
385 
386  if (fasterCode)
387  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
388 
389  // SI - For atom. deexc. tagging - 23/05/2017
390 
391  G4int Z = 8;
392 
393  G4ThreeVector deltaDirection =
394  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
395  Z, ionizationShell,
396  couple->GetMaterial());
397 
398  if (secondaryKinetic>0)
399  {
400  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
401  fvect->push_back(dp);
402  }
403 
404  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
405 
406  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
407  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
408  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
409  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
410  finalPx /= finalMomentum;
411  finalPy /= finalMomentum;
412  finalPz /= finalMomentum;
413 
414  G4ThreeVector direction;
415  direction.set(finalPx,finalPy,finalPz);
416 
418 
419  // AM: sample deexcitation
420  // here we assume that H_{2}O electronic levels are the same as Oxygen.
421  // this can be considered true with a rough 10% error in energy on K-shell,
422 
423  G4int secNumberInit = 0;// need to know at a certain point the energy of secondaries
424  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
425 
427  {
429 
430  if (ionizationShell <5 && ionizationShell >1)
431  {
432  as = G4AtomicShellEnumerator(4-ionizationShell);
433  }
434  else if (ionizationShell <2)
435  {
436  as = G4AtomicShellEnumerator(3);
437  }
438 
439  // FOR DEBUG ONLY
440  // if (ionizationShell == 4) {
441  //
442  // G4cout << "Z: " << Z << " as: " << as
443  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
444  // G4cout << "Press <Enter> key to continue..." << G4endl;
445  // G4cin.ignore();
446  // }
447 
448  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
449  secNumberInit = fvect->size();
450  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
451  secNumberFinal = fvect->size();
452  }
453 
454  // Note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
455 
456  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
457  G4double deexSecEnergy = 0;
458  for (G4int j=secNumberInit; j < secNumberFinal; j++)
459  {
460  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
461  }
462 
463  if (!statCode)
464  {
466  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
467  }
468  else
469  {
472  }
473 
474  // TEST //////////////////////////
475  // if (secondaryKinetic<0) abort();
476  // if (scatteredEnergy<0) abort();
477  // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
478  // if (k-scatteredEnergy<0) abort();
480 
481  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
483  ionizationShell,
484  theIncomingTrack);
485  }
486 
487 }
488 
489 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
490 
491 G4double
494  G4double k,
495  G4int shell)
496 {
497  // G4cout << "*** SLOW computation for "
498  // << " " << particleDefinition->GetParticleName() << G4endl;
499 
500  if(particleDefinition == G4Electron::ElectronDefinition())
501  {
502  G4double maximumEnergyTransfer = 0.;
503  if((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
504  maximumEnergyTransfer = k;
505  else
506  maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell))/ 2.;
507 
508  // SI : original method
509  /*
510  G4double crossSectionMaximum = 0.;
511  for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
512  {
513  G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
514  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
515  }
516  */
517 
518  // SI : alternative method
519  G4double crossSectionMaximum = 0.;
520 
521  G4double minEnergy = waterStructure.IonisationEnergy(shell);
522  G4double maxEnergy = maximumEnergyTransfer;
523  G4int nEnergySteps = 50;
524 
525  G4double value(minEnergy);
526  G4double stpEnergy(std::pow(maxEnergy / value,
527  1. / static_cast<G4double>(nEnergySteps - 1)));
528  G4int step(nEnergySteps);
529  while(step > 0)
530  {
531  step--;
532  G4double differentialCrossSection =
533  DifferentialCrossSection(particleDefinition,
534  k / eV,
535  value / eV,
536  shell);
537  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum =
538  differentialCrossSection;
539  value *= stpEnergy;
540  }
541  //
542 
543  G4double secondaryElectronKineticEnergy = 0.;
544  do
545  {
546  secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
547  }while(G4UniformRand()*crossSectionMaximum >
548  DifferentialCrossSection(particleDefinition, k/eV,
549  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
550 
551  return secondaryElectronKineticEnergy;
552 
553  }
554 
555  return 0;
556 }
557 
558 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
559 
560 // The following section is not used anymore but is kept for memory
561 // GetAngularDistribution()->SampleDirectionForShell is used instead
562 
563 /*
564  void G4DNAEmfietzoglouIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
565  G4double k,
566  G4double secKinetic,
567  G4double & cosTheta,
568  G4double & phi )
569  {
570  if (particleDefinition == G4Electron::ElectronDefinition())
571  {
572  phi = twopi * G4UniformRand();
573  if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
574  else if (secKinetic <= 200.*eV)
575  {
576  if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
577  else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
578  }
579  else
580  {
581  G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
582  cosTheta = std::sqrt(1.-sin2O);
583  }
584  }
585 
586  else if (particleDefinition == G4Proton::ProtonDefinition())
587  {
588  G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
589  phi = twopi * G4UniformRand();
590 
591  // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
592 
593  // Restriction below 100 eV from Emfietzoglou (2000)
594 
595  if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
596  else cosTheta = (2.*G4UniformRand())-1.;
597 
598  }
599  }
600  */
601 
602 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
604  G4double k,
605  G4double energyTransfer,
606  G4int ionizationLevelIndex)
607 {
608  G4double sigma = 0.;
609 
610  if(energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)/eV)
611  {
612  G4double valueT1 = 0;
613  G4double valueT2 = 0;
614  G4double valueE21 = 0;
615  G4double valueE22 = 0;
616  G4double valueE12 = 0;
617  G4double valueE11 = 0;
618 
619  G4double xs11 = 0;
620  G4double xs12 = 0;
621  G4double xs21 = 0;
622  G4double xs22 = 0;
623 
624  if(particleDefinition == G4Electron::ElectronDefinition())
625  {
626  // k should be in eV and energy transfer eV also
627 
628  std::vector<G4double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
629  eTdummyVec.end(),
630  k);
631 
632  std::vector<G4double>::iterator t1 = t2 - 1;
633 
634  // SI : the following condition avoids situations where energyTransfer >last vector element
635  // added strict limitations (09/08/2017)
636  if(energyTransfer < eVecm[(*t1)].back() &&
637  energyTransfer < eVecm[(*t2)].back())
638  {
639  std::vector<G4double>::iterator e12 =
640  std::upper_bound(eVecm[(*t1)].begin(),
641  eVecm[(*t1)].end(),
642  energyTransfer);
643  std::vector<G4double>::iterator e11 = e12 - 1;
644 
645  std::vector<G4double>::iterator e22 =
646  std::upper_bound(eVecm[(*t2)].begin(),
647  eVecm[(*t2)].end(),
648  energyTransfer);
649  std::vector<G4double>::iterator e21 = e22 - 1;
650 
651  valueT1 = *t1;
652  valueT2 = *t2;
653  valueE21 = *e21;
654  valueE22 = *e22;
655  valueE12 = *e12;
656  valueE11 = *e11;
657 
658  xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
659  xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
660  xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
661  xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
662 
663  //G4cout << "-------------------" << G4endl;
664  //G4cout << "ionizationLevelIndex=" << ionizationLevelIndex << G4endl;
665  //G4cout << "valueT1/eV=" << valueT1 << " valueT2/eV=" << valueT2 << G4endl;
666  //G4cout << "valueE11/eV=" << valueE11 << " valueE12/eV=" << valueE12 << " valueE21/eV=" << valueE21 << " valueE22/eV=" << valueE22 << G4endl;
667  //G4cout << "xs11=" << xs11 / ((1.e-22 / 3.343) * m*m) << G4endl;
668  //G4cout << "xs12=" << xs12 / ((1.e-22 / 3.343) * m*m) << G4endl;
669  //G4cout << "xs21=" << xs21 / ((1.e-22 / 3.343) * m*m) << G4endl;
670  //G4cout << "xs22=" << xs22 / ((1.e-22 / 3.343) * m*m) << G4endl;
671  //G4cout << "###################" << G4endl;
672 
673  }
674 
675  }
676 
677  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
678  if(xsProduct != 0.)
679  {
680  sigma = QuadInterpolator(valueE11,
681  valueE12,
682  valueE21,
683  valueE22,
684  xs11,
685  xs12,
686  xs21,
687  xs22,
688  valueT1,
689  valueT2,
690  k,
691  energyTransfer);
692  }
693 
694  }
695 
696  return sigma;
697 }
698 
699 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
700 
702  G4double e2,
703  G4double e,
704  G4double xs1,
705  G4double xs2)
706 {
707 
708  G4double value = 0.;
709 
710  // Log-log interpolation by default
711 
712  if(e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
713  && !fasterCode)
714  {
715  G4double a = (std::log10(xs2) - std::log10(xs1))
716  / (std::log10(e2) - std::log10(e1));
717  G4double b = std::log10(xs2) - a * std::log10(e2);
718  G4double sigma = a * std::log10(e) + b;
719  value = (std::pow(10., sigma));
720  }
721 
722  // Switch to lin-lin interpolation
723  /*
724  if ((e2-e1)!=0)
725  {
726  G4double d1 = xs1;
727  G4double d2 = xs2;
728  value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
729  }
730  */
731 
732  // Switch to log-lin interpolation for faster code
733  if((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
734  {
735  G4double d1 = std::log10(xs1);
736  G4double d2 = std::log10(xs2);
737  value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
738  }
739 
740  // Switch to lin-lin interpolation for faster code
741  // in case one of xs1 or xs2 (=cum proba) value is zero
742 
743  if((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
744  {
745  G4double d1 = xs1;
746  G4double d2 = xs2;
747  value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
748  }
749 
750  /*
751  G4cout
752  << e1 << " "
753  << e2 << " "
754  << e << " "
755  << xs1 << " "
756  << xs2 << " "
757  << value
758  << G4endl;
759  */
760 
761  return value;
762 }
763 
764 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
765 
767  G4double e12,
768  G4double e21,
769  G4double e22,
770  G4double xs11,
771  G4double xs12,
772  G4double xs21,
773  G4double xs22,
774  G4double t1,
775  G4double t2,
776  G4double t,
777  G4double e)
778 {
779  G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
780  G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
782  t2,
783  t,
784  interpolatedvalue1,
785  interpolatedvalue2);
786 
787  return value;
788 }
789 
790 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
791 
793  const G4String& particle)
794 {
795  G4int level = 0;
796 
797  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
798  pos = tableData.find(particle);
799 
800  if(pos != tableData.end())
801  {
802  G4DNACrossSectionDataSet* table = pos->second;
803 
804  if(table != 0)
805  {
806  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
807  const size_t n(table->NumberOfComponents());
808  size_t i(n);
809  G4double value = 0.;
810 
811  while(i > 0)
812  {
813  i--;
814  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
815  value += valuesBuffer[i];
816  }
817 
818  value *= G4UniformRand();
819 
820  i = n;
821 
822  while(i > 0)
823  {
824  i--;
825 
826  if(valuesBuffer[i] > value)
827  {
828  delete[] valuesBuffer;
829  return i;
830  }
831  value -= valuesBuffer[i];
832  }
833 
834  if(valuesBuffer) delete[] valuesBuffer;
835 
836  }
837  }
838  else
839  {
840  G4Exception("G4DNAEmfietzoglouIonisationModel::RandomSelect",
841  "em0002",
843  "Model not applicable to particle type.");
844  }
845 
846  return level;
847 }
848 
849 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
850 
852  G4double k,
853  G4int shell)
854 {
855  //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
856 
857  G4double secondaryElectronKineticEnergy = 0.;
858 
859  secondaryElectronKineticEnergy = RandomTransferedEnergy(particleDefinition,
860  k / eV,
861  shell)
862  * eV
864 
865  //G4cout << RandomTransferedEnergy(particleDefinition, k/eV, shell) << G4endl;
866  if(secondaryElectronKineticEnergy < 0.) return 0.;
867 
868  return secondaryElectronKineticEnergy;
869 }
870 
871 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
872 
874  G4double k,
875  G4int ionizationLevelIndex)
876 {
877 
878  G4double random = G4UniformRand();
879 
880  G4double nrj = 0.;
881 
882  G4double valueK1 = 0;
883  G4double valueK2 = 0;
884  G4double valuePROB21 = 0;
885  G4double valuePROB22 = 0;
886  G4double valuePROB12 = 0;
887  G4double valuePROB11 = 0;
888 
889  G4double nrjTransf11 = 0;
890  G4double nrjTransf12 = 0;
891  G4double nrjTransf21 = 0;
892  G4double nrjTransf22 = 0;
893 
894  if (particleDefinition == G4Electron::ElectronDefinition())
895  {
896  // k should be in eV
897  std::vector<G4double>::iterator k2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
898 
899  std::vector<G4double>::iterator k1 = k2-1;
900 
901  /*
902  G4cout << "----> k=" << k
903  << " " << *k1
904  << " " << *k2
905  << " " << random
906  << " " << ionizationLevelIndex
907  << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
908  << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
909  << G4endl;
910  */
911 
912  // SI : the following condition avoids situations where random >last vector element
913  if ( random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
914  && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back() )
915 
916  {
917  std::vector<G4double>::iterator prob12 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
918  eProbaShellMap[ionizationLevelIndex][(*k1)].end(), random);
919 
920  std::vector<G4double>::iterator prob11 = prob12-1;
921 
922  std::vector<G4double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
923  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
924 
925  std::vector<G4double>::iterator prob21 = prob22-1;
926 
927  valueK1 =*k1;
928  valueK2 =*k2;
929  valuePROB21 =*prob21;
930  valuePROB22 =*prob22;
931  valuePROB12 =*prob12;
932  valuePROB11 =*prob11;
933 
934  /*
935  G4cout << " " << random << " " << valuePROB11 << " "
936  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
937  */
938 
939  nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
940  nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
941  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
942  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
943 
944  /*
945  G4cout << " " << ionizationLevelIndex << " "
946  << random << " " <<valueK1 << " " << valueK2 << G4endl;
947 
948  G4cout << " " << random << " " << nrjTransf11 << " "
949  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
950  */
951 
952  }
953 
954  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
955 
956  if ( random > eProbaShellMap[ionizationLevelIndex][(*k1)].back() )
957 
958  {
959  std::vector<G4double>::iterator prob22 = std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
960  eProbaShellMap[ionizationLevelIndex][(*k2)].end(), random);
961 
962  std::vector<G4double>::iterator prob21 = prob22-1;
963 
964  valueK1 =*k1;
965  valueK2 =*k2;
966  valuePROB21 =*prob21;
967  valuePROB22 =*prob22;
968 
969  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
970 
971  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
972  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
973 
974  G4double interpolatedvalue2 = Interpolate(valuePROB21, valuePROB22, random, nrjTransf21, nrjTransf22);
975 
976  // zeros are explicitely set
977 
978  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
979 
980  /*
981  G4cout << " " << ionizationLevelIndex << " "
982  << random << " " <<valueK1 << " " << valueK2 << G4endl;
983 
984  G4cout << " " << random << " " << nrjTransf11 << " "
985  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
986 
987  G4cout << "ici" << " " << value << G4endl;
988  */
989 
990  return value;
991  }
992 
993  }
994 
995  // End electron
996 
997  G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21 * nrjTransf22;
998 
999  //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
1000 
1001  if (nrjTransfProduct != 0.)
1002  {
1003  nrj = QuadInterpolator( valuePROB11, valuePROB12,
1004  valuePROB21, valuePROB22,
1005  nrjTransf11, nrjTransf12,
1006  nrjTransf21, nrjTransf22,
1007  valueK1, valueK2,
1008  k, random);
1009  }
1010 
1011  //G4cout << nrj << endl;
1012 
1013  return nrj;
1014 }
void set(double x, double y, double z)
virtual G4bool LoadData(const G4String &argFileName)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static const G4double pos
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
TTree * t1
Definition: plottest35.C:26
G4VAtomDeexcitation * AtomDeexcitation()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
virtual size_t NumberOfComponents(void) const
static constexpr double keV
Definition: G4SIunits.hh:216
size_t GetIndex() const
Definition: G4Material.hh:262
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const G4String & GetParticleName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4int RandomSelect(G4double energy, const G4String &particle)
double z() const
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)
G4AtomicShellEnumerator
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:608
G4double GetPDGMass() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4double RandomTransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
Float_t Z
static constexpr double m
Definition: G4SIunits.hh:129
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:763
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetDefinition() const
static const G4double d2
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static constexpr double electron_mass_c2
static G4DNAMolecularMaterial * Instance()
static const G4double d1
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:53
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4DNAEmfietzoglouWaterIonisationStructure waterStructure
static constexpr double eV
Definition: G4SIunits.hh:215
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:582
const G4Track * GetCurrentTrack() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
Hep3Vector unit() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
TTree * t2
Definition: plottest35.C:36
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
int G4int
Definition: G4Types.hh:78
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
const std::vector< G4double > * fpMolWaterDensity
G4double GetKineticEnergy() const
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4DNAEmfietzoglouIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
double x() const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
static G4LossTableManager * Instance()
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:575
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
Char_t n[5]
const G4Material * GetMaterial() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
double y() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
static G4DNAChemistryManager * Instance()
virtual G4double FindValue(G4double x, G4int componentId=0) const =0