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