Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermorePhotoElectricModel.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: G4LivermorePhotoElectricModel.cc 79186 2014-02-20 09:20:02Z gcosmo $
27 //
28 //
29 // Author: Sebastien Incerti
30 // 30 October 2008
31 // on base of G4LowEnergyPhotoElectric developed by A.Forti and M.G.Pia
32 //
33 // 22 Oct 2012 A & V Ivanchenko Migration data structure to G4PhysicsVector
34 //
35 
37 #include "G4SystemOfUnits.hh"
38 #include "G4LossTableManager.hh"
39 #include "G4Electron.hh"
40 #include "G4Gamma.hh"
42 #include "G4CrossSectionHandler.hh"
43 #include "G4LPhysicsFreeVector.hh"
44 #include "G4VAtomDeexcitation.hh"
46 #include "G4AtomicShell.hh"
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
50 G4LPhysicsFreeVector* G4LivermorePhotoElectricModel::fCrossSection[] = {0};
51 G4LPhysicsFreeVector* G4LivermorePhotoElectricModel::fCrossSectionLE[] = {0};
52 std::vector<G4double>* G4LivermorePhotoElectricModel::fParam[] = {0};
53 G4int G4LivermorePhotoElectricModel::fNShells[] = {0};
54 G4int G4LivermorePhotoElectricModel::fNShellsUsed[] = {0};
55 G4ElementData* G4LivermorePhotoElectricModel::fShellCrossSection = 0;
56 
57 using namespace std;
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60 
62  const G4String& nam)
63  : G4VEmModel(nam),fParticleChange(0),maxZ(99),
64  nShellLimit(100),fDeexcitationActive(false),isInitialised(false),
65  fAtomDeexcitation(0)
66 {
67  verboseLevel= 0;
68  // Verbosity scale:
69  // 0 = nothing
70  // 1 = warning for energy non-conservation
71  // 2 = details of energy budget
72  // 3 = calculation of cross sections, file openings, sampling of atoms
73  // 4 = entering in methods
74 
75  theGamma = G4Gamma::Gamma();
76  theElectron = G4Electron::Electron();
77 
78  // default generator
80 
81  if(verboseLevel>0) {
82  G4cout << "Livermore PhotoElectric is constructed "
83  << " nShellLimit= " << nShellLimit << G4endl;
84  }
85 
86  //Mark this model as "applicable" for atomic deexcitation
87  SetDeexcitationFlag(true);
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
94  if(IsMaster()) {
95  delete fShellCrossSection;
96  for(G4int i=0; i<maxZ; ++i) {
97  delete fParam[i];
98  }
99  }
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
104 void
106  const G4DataVector&)
107 {
108  if (verboseLevel > 2) {
109  G4cout << "Calling G4LivermorePhotoElectricModel::Initialise()" << G4endl;
110  }
111 
112  if(IsMaster()) {
113 
114  if(!fShellCrossSection) { fShellCrossSection = new G4ElementData(); }
115 
116  char* path = getenv("G4LEDATA");
117 
118  G4ProductionCutsTable* theCoupleTable =
120  G4int numOfCouples = theCoupleTable->GetTableSize();
121 
122  for(G4int i=0; i<numOfCouples; ++i) {
123  const G4MaterialCutsCouple* couple =
124  theCoupleTable->GetMaterialCutsCouple(i);
125  const G4Material* material = couple->GetMaterial();
126  const G4ElementVector* theElementVector = material->GetElementVector();
127  G4int nelm = material->GetNumberOfElements();
128 
129  for (G4int j=0; j<nelm; ++j) {
130  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
131  if(Z < 1) { Z = 1; }
132  else if(Z > maxZ) { Z = maxZ; }
133  if(!fCrossSection[Z]) { ReadData(Z, path); }
134  }
135  }
136  }
137  //
138  if (verboseLevel > 2) {
139  G4cout << "Loaded cross section files for LivermorePhotoElectric model"
140  << G4endl;
141  }
142  if(!isInitialised) {
143  isInitialised = true;
145 
146  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
147  }
148  fDeexcitationActive = false;
149  if(fAtomDeexcitation) {
150  fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
151  }
152 
153  if (verboseLevel > 0) {
154  G4cout << "LivermorePhotoElectric model is initialized " << G4endl
155  << G4endl;
156  }
157 }
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160 
162  const G4ParticleDefinition*,
164  G4double ZZ, G4double,
166 {
167  if (verboseLevel > 3) {
168  G4cout << "G4LivermorePhotoElectricModel::ComputeCrossSectionPerAtom():"
169  << " Z= " << ZZ << " R(keV)= " << energy/keV << G4endl;
170  }
171  G4double cs = 0.0;
172  G4double gammaEnergy = energy;
173 
174  G4int Z = G4lrint(ZZ);
175  if(Z < 1 || Z >= maxZ) { return cs; }
176 
177  // if element was not initialised
178  // do initialisation safely for MT mode
179  if(!fCrossSection[Z]) {
180  InitialiseForElement(0, Z);
181  if(!fCrossSection[Z]) { return cs; }
182  }
183 
184  G4int idx = fNShells[Z]*6 - 4;
185  if (gammaEnergy <= (*(fParam[Z]))[idx-1]) { return cs; }
186 
187  G4double x1 = 1.0/gammaEnergy;
188  G4double x2 = x1*x1;
189  G4double x3 = x2*x1;
190 
191  // parameterisation
192  if(gammaEnergy >= (*(fParam[Z]))[0]) {
193  G4double x4 = x2*x2;
194  cs = x1*((*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
195  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
196  + x4*(*(fParam[Z]))[idx+4]);
197  // high energy part
198  } else if(gammaEnergy >= (*(fParam[Z]))[1]) {
199  cs = x3*(fCrossSection[Z])->Value(gammaEnergy);
200 
201  // low energy part
202  } else {
203  cs = x3*(fCrossSectionLE[Z])->Value(gammaEnergy);
204  }
205  if (verboseLevel > 1) {
206  G4cout << "LivermorePhotoElectricModel: E(keV)= " << gammaEnergy/keV
207  << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
208  }
209  return cs;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213 
214 void
216  std::vector<G4DynamicParticle*>* fvect,
217  const G4MaterialCutsCouple* couple,
218  const G4DynamicParticle* aDynamicGamma,
219  G4double,
220  G4double)
221 {
222  G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
223  if (verboseLevel > 3) {
224  G4cout << "G4LivermorePhotoElectricModel::SampleSecondaries() Egamma(keV)= "
225  << gammaEnergy/keV << G4endl;
226  }
227 
228  // kill incident photon
231 
232  // Returns the normalized direction of the momentum
233  G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
234 
235  // Select randomly one element in the current material
236  //G4cout << "Select random atom Egamma(keV)= " << gammaEnergy/keV << G4endl;
237  const G4Element* elm = SelectRandomAtom(couple->GetMaterial(),theGamma,
238  gammaEnergy);
239  G4int Z = G4lrint(elm->GetZ());
240 
241  // Select the ionised shell in the current atom according to shell
242  // cross sections
243  // G4cout << "Select random shell Z= " << Z << G4endl;
244 
245  if(Z >= maxZ) { Z = maxZ-1; }
246 
247  // element was not initialised gamma should be absorbed
248  if(!fCrossSection[Z]) {
250  return;
251  }
252 
253  // shell index
254  size_t shellIdx = 0;
255  size_t nn = fNShellsUsed[Z];
256 
257  if(nn > 1) {
258  if(gammaEnergy >= (*(fParam[Z]))[0]) {
259  G4double x1 = 1.0/gammaEnergy;
260  G4double x2 = x1*x1;
261  G4double x3 = x2*x1;
262  G4double x4 = x3*x1;
263  G4int idx = nn*6 - 4;
264  // when do sampling common factors are not taken into account
265  // so cross section is not real
266  G4double cs0 = G4UniformRand()*((*(fParam[Z]))[idx]
267  + x1*(*(fParam[Z]))[idx+1]
268  + x2*(*(fParam[Z]))[idx+2]
269  + x3*(*(fParam[Z]))[idx+3]
270  + x4*(*(fParam[Z]))[idx+4]);
271  for(shellIdx=0; shellIdx<nn; ++shellIdx) {
272  idx = shellIdx*6 + 2;
273  if(gammaEnergy > (*(fParam[Z]))[idx-1]) {
274  G4double cs = (*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
275  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
276  + x4*(*(fParam[Z]))[idx+4];
277  if(cs >= cs0) { break; }
278  }
279  }
280  if(shellIdx >= nn) { shellIdx = nn-1; }
281 
282  } else {
283 
284  // when do sampling common factors are not taken into account
285  // so cross section is not real
286  G4double cs = G4UniformRand();
287 
288  if(gammaEnergy >= (*(fParam[Z]))[1]) {
289  cs *= (fCrossSection[Z])->Value(gammaEnergy);
290  } else {
291  cs *= (fCrossSectionLE[Z])->Value(gammaEnergy);
292  }
293 
294  for(size_t j=0; j<nn; ++j) {
295  shellIdx = (size_t)fShellCrossSection->GetComponentID(Z, j);
296  if(gammaEnergy > (*(fParam[Z]))[6*shellIdx+1]) {
297  cs -= fShellCrossSection->GetValueForComponent(Z, j, gammaEnergy);
298  }
299  if(cs <= 0.0 || j+1 == nn) { break; }
300  }
301  }
302  }
303 
304  G4double bindingEnergy = (*(fParam[Z]))[shellIdx*6 + 1];
305  //G4cout << "Z= " << Z << " shellIdx= " << shellIdx
306  // << " nShells= " << fNShells[Z]
307  // << " Ebind(keV)= " << bindingEnergy/keV
308  // << " Egamma(keV)= " << gammaEnergy/keV << G4endl;
309 
310  const G4AtomicShell* shell = 0;
311 
312  // no de-excitation from the last shell
313  if(fDeexcitationActive && shellIdx + 1 < nn) {
315  shell = fAtomDeexcitation->GetAtomicShell(Z, as);
316  }
317 
318  // If binding energy of the selected shell is larger than photon energy
319  // do not generate secondaries
320  if(gammaEnergy < bindingEnergy) {
322  return;
323  }
324 
325  // Primary outcoming electron
326  G4double eKineticEnergy = gammaEnergy - bindingEnergy;
327  G4double edep = bindingEnergy;
328 
329  // Calculate direction of the photoelectron
330  G4ThreeVector electronDirection =
331  GetAngularDistribution()->SampleDirection(aDynamicGamma,
332  eKineticEnergy,
333  shellIdx,
334  couple->GetMaterial());
335 
336  // The electron is created
337  G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
338  electronDirection,
339  eKineticEnergy);
340  fvect->push_back(electron);
341 
342  // Sample deexcitation
343  if(shell) {
344  G4int index = couple->GetIndex();
345  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
346  size_t nbefore = fvect->size();
347 
348  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
349  size_t nafter = fvect->size();
350  if(nafter > nbefore) {
351  for (size_t j=nbefore; j<nafter; ++j) {
352  edep -= ((*fvect)[j])->GetKineticEnergy();
353  }
354  }
355  }
356  }
357  // energy balance - excitation energy left
358  if(edep > 0.0) {
360  }
361 }
362 
363 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
364 
365 void
366 G4LivermorePhotoElectricModel::ReadData(G4int Z, const char* path)
367 {
368  if (verboseLevel > 1)
369  {
370  G4cout << "Calling ReadData() of G4LivermoreGammaConversionModel"
371  << G4endl;
372  }
373 
374  if(fCrossSection[Z]) { return; }
375 
376  const char* datadir = path;
377 
378  if(!datadir)
379  {
380  datadir = getenv("G4LEDATA");
381  if(!datadir)
382  {
383  G4Exception("G4LivermorePhotoElectricModel::ReadData()",
384  "em0006",FatalException,
385  "Environment variable G4LEDATA not defined");
386  return;
387  }
388  }
389 
390  // spline for photoeffect total x-section above K-shell
391  fCrossSection[Z] = new G4LPhysicsFreeVector();
392  fCrossSection[Z]->SetSpline(true);
393 
394  std::ostringstream ost;
395  ost << datadir << "/livermore/phot/pe-cs-" << Z <<".dat";
396  std::ifstream fin(ost.str().c_str());
397  if( !fin.is_open()) {
399  ed << "G4LivermorePhotoElectricModel data file <" << ost.str().c_str()
400  << "> is not opened!" << G4endl;
401  G4Exception("G4LivermorePhotoElectricModel::ReadData()",
402  "em0003",FatalException,
403  ed,"G4LEDATA version should be G4EMLOW6.32 or later.");
404  return;
405  } else {
406  if(verboseLevel > 3) { G4cout << "File " << ost.str().c_str()
407  << " is opened by G4LivermorePhotoElectricModel" << G4endl;}
408  fCrossSection[Z]->Retrieve(fin, true);
409  fCrossSection[Z]->ScaleVector(MeV, barn);
410  fin.close();
411  }
412 
413  fParam[Z] = new std::vector<G4double>;
414 
415  // read fit parameters
416  G4int n1 = 0;
417  G4int n2 = 0;
418  G4double x;
419  std::ostringstream ost1;
420  ost1 << datadir << "/livermore/phot/pe-" << Z <<".dat";
421  std::ifstream fin1(ost1.str().c_str());
422  if( !fin1.is_open()) {
424  ed << "G4LivermorePhotoElectricModel data file <" << ost1.str().c_str()
425  << "> is not opened!" << G4endl;
426  G4Exception("G4LivermorePhotoElectricModel::ReadData()",
427  "em0003",FatalException,
428  ed,"G4LEDATA version should be G4EMLOW6.32 or later.");
429  return;
430  } else {
431  if(verboseLevel > 3) {
432  G4cout << "File " << ost1.str().c_str()
433  << " is opened by G4LivermorePhotoElectricModel" << G4endl;
434  }
435  fin1 >> n1;
436  if(fin1.fail()) { return; }
437  if(0 > n1 || n1 >= INT_MAX) { n1 = 0; }
438 
439  fin1 >> n2;
440  if(fin1.fail()) { return; }
441  if(0 > n2 || n2 >= INT_MAX) { n2 = 0; }
442 
443  fin1 >> x;
444  if(fin1.fail()) { return; }
445 
446  fNShells[Z] = n1;
447  fParam[Z]->reserve(6*n1+1);
448  fParam[Z]->push_back(x*MeV);
449  for(G4int i=0; i<n1; ++i) {
450  for(G4int j=0; j<6; ++j) {
451  fin1 >> x;
452  if(0 == j) { x *= MeV; }
453  else { x *= barn; }
454  fParam[Z]->push_back(x);
455  }
456  }
457  fin1.close();
458  }
459  // there is a possibility to used only main shells
460  if(nShellLimit < n2) { n2 = nShellLimit; }
461  fShellCrossSection->InitialiseForComponent(Z, n2);
462  fNShellsUsed[Z] = n2;
463 
464  if(1 < n2) {
465  std::ostringstream ost2;
466  ost2 << datadir << "/livermore/phot/pe-ss-cs-" << Z <<".dat";
467  std::ifstream fin2(ost2.str().c_str());
468  if( !fin2.is_open()) {
470  ed << "G4LivermorePhotoElectricModel data file <" << ost2.str().c_str()
471  << "> is not opened!" << G4endl;
472  G4Exception("G4LivermorePhotoElectricModel::ReadData()",
473  "em0003",FatalException,
474  ed,"G4LEDATA version should be G4EMLOW6.32 or later.");
475  return;
476  } else {
477  if(verboseLevel > 3) {
478  G4cout << "File " << ost2.str().c_str()
479  << " is opened by G4LivermorePhotoElectricModel" << G4endl;
480  }
481 
482  G4int n3, n4;
483  G4double y;
484  for(G4int i=0; i<n2; ++i) {
485  fin2 >> x >> y >> n3 >> n4;
487  for(G4int j=0; j<n3; ++j) {
488  fin2 >> x >> y;
489  v->PutValues(j, x*MeV, y*barn);
490  }
491  fShellCrossSection->AddComponent(Z, n4, v);
492  }
493  fin2.close();
494  }
495  }
496 
497  // no spline for photoeffect total x-section below K-shell
498  if(1 < fNShells[Z]) {
499  fCrossSectionLE[Z] = new G4LPhysicsFreeVector();
500 
501  std::ostringstream ost3;
502  ost3 << datadir << "/livermore/phot/pe-le-cs-" << Z <<".dat";
503  std::ifstream fin3(ost3.str().c_str());
504  if( !fin3.is_open()) {
506  ed << "G4LivermorePhotoElectricModel data file <" << ost3.str().c_str()
507  << "> is not opened!" << G4endl;
508  G4Exception("G4LivermorePhotoElectricModel::ReadData()",
509  "em0003",FatalException,
510  ed,"G4LEDATA version should be G4EMLOW6.32 or later.");
511  return;
512  } else {
513  if(verboseLevel > 3) {
514  G4cout << "File " << ost3.str().c_str()
515  << " is opened by G4LivermorePhotoElectricModel" << G4endl;
516  }
517  fCrossSectionLE[Z]->Retrieve(fin3, true);
518  fCrossSectionLE[Z]->ScaleVector(MeV, barn);
519  fin3.close();
520  }
521  }
522 }
523 
524 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
525 
526 #include "G4AutoLock.hh"
527 namespace { G4Mutex LivermorePhotoElectricModelMutex = G4MUTEX_INITIALIZER; }
528 
530  const G4ParticleDefinition*, G4int Z)
531 {
532  G4AutoLock l(&LivermorePhotoElectricModelMutex);
533  // G4cout << "G4LivermorePhotoElectricModel::InitialiseForElement Z= "
534  // << Z << G4endl;
535  if(!fCrossSection[Z]) { ReadData(Z); }
536  l.unlock();
537 }
538 
539 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void PutValues(size_t binNumber, G4double binValue, G4double dataValue)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
std::vector< G4Element * > G4ElementVector
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4LivermorePhotoElectricModel(const G4String &nam="LivermorePhElectric")
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetSpline(G4bool)
string material
Definition: eplot.py:19
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4int GetComponentID(G4int Z, size_t idx)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
const G4ThreeVector & GetMomentumDirection() const
virtual void ScaleVector(G4double factorE, G4double factorV)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define INT_MAX
Definition: templates.hh:111
static G4ProductionCutsTable * GetProductionCutsTable()
G4int G4Mutex
Definition: G4Threading.hh:156
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double GetValueForComponent(G4int Z, size_t idx, G4double kinEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double bindingEnergy(G4int A, G4int Z)
void InitialiseForComponent(G4int Z, G4int nComponents=0)
G4ParticleChangeForGamma * fParticleChange
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
const G4Material * GetMaterial() const