Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes
G4LivermorePhotoElectricModel Class Reference

#include <G4LivermorePhotoElectricModel.hh>

Inheritance diagram for G4LivermorePhotoElectricModel:
G4VEmModel

Public Member Functions

 G4LivermorePhotoElectricModel (const G4String &nam="LivermorePhElectric")
 
virtual ~G4LivermorePhotoElectricModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
void SetLimitNumberOfShells (G4int)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 50 of file G4LivermorePhotoElectricModel.hh.

Constructor & Destructor Documentation

G4LivermorePhotoElectricModel::G4LivermorePhotoElectricModel ( const G4String nam = "LivermorePhElectric")

Definition at line 61 of file G4LivermorePhotoElectricModel.cc.

References G4Electron::Electron(), G4cout, G4endl, G4Gamma::Gamma(), G4VEmModel::SetAngularDistribution(), and G4VEmModel::SetDeexcitationFlag().

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 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4ParticleChangeForGamma * fParticleChange
G4LivermorePhotoElectricModel::~G4LivermorePhotoElectricModel ( )
virtual

Definition at line 92 of file G4LivermorePhotoElectricModel.cc.

References G4VEmModel::IsMaster().

93 {
94  if(IsMaster()) {
95  delete fShellCrossSection;
96  for(G4int i=0; i<maxZ; ++i) {
97  delete fParam[i];
98  }
99  }
100 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
int G4int
Definition: G4Types.hh:78

Member Function Documentation

G4double G4LivermorePhotoElectricModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 161 of file G4LivermorePhotoElectricModel.cc.

References python.hepunit::barn, energy(), G4cout, G4endl, G4lrint(), InitialiseForElement(), python.hepunit::keV, and G4VEmModel::Value().

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 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
int G4int
Definition: G4Types.hh:78
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
void G4LivermorePhotoElectricModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 105 of file G4LivermorePhotoElectricModel.cc.

References G4LossTableManager::AtomDeexcitation(), fParticleChange, G4cout, G4endl, G4Material::GetElementVector(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4VEmModel::GetParticleChangeForGamma(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4LossTableManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), G4VEmModel::IsMaster(), and eplot::material.

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 }
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::vector< G4Element * > G4ElementVector
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * fParticleChange
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const
void G4LivermorePhotoElectricModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 529 of file G4LivermorePhotoElectricModel.cc.

References G4TemplateAutoLock< M, L, U >::unlock().

Referenced by ComputeCrossSectionPerAtom().

531 {
532  G4AutoLock l(&LivermorePhotoElectricModelMutex);
533  // G4cout << "G4LivermorePhotoElectricModel::InitialiseForElement Z= "
534  // << Z << G4endl;
535  if(!fCrossSection[Z]) { ReadData(Z); }
536  l.unlock();
537 }
void G4LivermorePhotoElectricModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 215 of file G4LivermorePhotoElectricModel.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4VAtomDeexcitation::CheckDeexcitationActiveRegion(), G4InuclParticleNames::electron, fParticleChange, fStopAndKill, G4cout, G4endl, G4lrint(), G4UniformRand, G4VAtomDeexcitation::GenerateParticles(), G4VEmModel::GetAngularDistribution(), G4VAtomDeexcitation::GetAtomicShell(), G4ElementData::GetComponentID(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4ElementData::GetValueForComponent(), G4Element::GetZ(), python.hepunit::keV, G4InuclParticleNames::nn, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and G4VEmModel::Value().

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 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4int GetComponentID(G4int Z, size_t idx)
#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
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetValueForComponent(G4int Z, size_t idx, G4double kinEnergy)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4double bindingEnergy(G4int A, G4int Z)
G4ParticleChangeForGamma * fParticleChange
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
const G4Material * GetMaterial() const
void G4LivermorePhotoElectricModel::SetLimitNumberOfShells ( G4int  n)
inline

Definition at line 112 of file G4LivermorePhotoElectricModel.hh.

References n.

113 {
114  nShellLimit = n;
115 }
const G4int n

Field Documentation

G4ParticleChangeForGamma* G4LivermorePhotoElectricModel::fParticleChange
protected

Definition at line 82 of file G4LivermorePhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


The documentation for this class was generated from the following files: