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