Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LivermorePhotoElectricModel Class Reference

#include <G4LivermorePhotoElectricModel.hh>

Inheritance diagram for G4LivermorePhotoElectricModel:
Collaboration diagram for G4LivermorePhotoElectricModel:

Public Member Functions

 G4LivermorePhotoElectricModel (const G4String &nam="LivermorePhElectric")
 
virtual ~G4LivermorePhotoElectricModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double energy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double energy, 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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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=nullptr)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 50 of file G4LivermorePhotoElectricModel.hh.

Constructor & Destructor Documentation

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

Definition at line 64 of file G4LivermorePhotoElectricModel.cc.

66  : G4VEmModel(nam),fParticleChange(nullptr),maxZ(99),
67  nShellLimit(100),fDeexcitationActive(false),isInitialised(false),
68  fAtomDeexcitation(nullptr)
69 {
70  verboseLevel= 0;
71  // Verbosity scale:
72  // 0 = nothing
73  // 1 = warning for energy non-conservation
74  // 2 = details of energy budget
75  // 3 = calculation of cross sections, file openings, sampling of atoms
76  // 4 = entering in methods
77 
78  theGamma = G4Gamma::Gamma();
79  theElectron = G4Electron::Electron();
80 
81  // default generator
83 
84  if(verboseLevel>0) {
85  G4cout << "Livermore PhotoElectric is constructed "
86  << " nShellLimit= " << nShellLimit << G4endl;
87  }
88 
89  //Mark this model as "applicable" for atomic deexcitation
90  SetDeexcitationFlag(true);
91  fSandiaCof.resize(4,0.0);
92  fCurrSection = 0.0;
93 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
G4ParticleChangeForGamma * fParticleChange

Here is the call graph for this function:

G4LivermorePhotoElectricModel::~G4LivermorePhotoElectricModel ( )
virtual

Definition at line 97 of file G4LivermorePhotoElectricModel.cc.

98 {
99  if(IsMaster()) {
100  delete fShellCrossSection;
101  for(G4int i=0; i<maxZ; ++i) {
102  delete fParam[i];
103  fParam[i] = 0;
104  delete fCrossSection[i];
105  fCrossSection[i] = 0;
106  delete fCrossSectionLE[i];
107  fCrossSectionLE[i] = 0;
108  }
109  }
110 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 205 of file G4LivermorePhotoElectricModel.cc.

210 {
211  if (verboseLevel > 3) {
212  G4cout << "G4LivermorePhotoElectricModel::ComputeCrossSectionPerAtom():"
213  << " Z= " << ZZ << " R(keV)= " << energy/keV << G4endl;
214  }
215  G4double cs = 0.0;
216  G4int Z = G4lrint(ZZ);
217  if(Z < 1 || Z >= maxZ) { return cs; }
218 
219  // if element was not initialised
220  // do initialisation safely for MT mode
221  if(!fCrossSection[Z]) {
222  InitialiseForElement(0, Z);
223  if(!fCrossSection[Z]) { return cs; }
224  }
225 
226  G4int idx = fNShells[Z]*6 - 4;
227  if (energy < (*(fParam[Z]))[idx-1]) { energy = (*(fParam[Z]))[idx-1]; }
228 
229  G4double x1 = 1.0/energy;
230  G4double x2 = x1*x1;
231  G4double x3 = x2*x1;
232 
233  // parameterisation
234  if(energy >= (*(fParam[Z]))[0]) {
235  G4double x4 = x2*x2;
236  cs = x1*((*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
237  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
238  + x4*(*(fParam[Z]))[idx+4]);
239  // high energy part
240  } else if(energy >= (*(fParam[Z]))[1]) {
241  cs = x3*(fCrossSection[Z])->Value(energy);
242 
243  // low energy part
244  } else {
245  cs = x3*(fCrossSectionLE[Z])->Value(energy);
246  }
247  if (verboseLevel > 1) {
248  G4cout << "LivermorePhotoElectricModel: E(keV)= " << energy/keV
249  << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
250  }
251  return cs;
252 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
G4double energy(const ThreeVector &p, const G4double m)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double barn
Definition: G4SIunits.hh:105
static constexpr double keV
Definition: G4SIunits.hh:216
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377

Here is the call graph for this function:

G4double G4LivermorePhotoElectricModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  energy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 176 of file G4LivermorePhotoElectricModel.cc.

181 {
182  fCurrSection = 0.0;
183  if(fWater && (material == fWater ||
184  material->GetBaseMaterial() == fWater)) {
185  if(energy <= fWaterEnergyLimit) {
186  fWater->GetSandiaTable()->GetSandiaCofWater(energy, fSandiaCof);
187 
188  G4double energy2 = energy*energy;
189  G4double energy3 = energy*energy2;
190  G4double energy4 = energy2*energy2;
191 
192  fCurrSection = material->GetDensity()*
193  (fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
194  fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4);
195  }
196  }
197  if(0.0 == fCurrSection) {
198  fCurrSection = G4VEmModel::CrossSectionPerVolume(material, p, energy);
199  }
200  return fCurrSection;
201 }
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:257
G4double GetDensity() const
Definition: G4Material.hh:180
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:229
G4double energy(const ThreeVector &p, const G4double m)
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
void GetSandiaCofWater(G4double energy, std::vector< G4double > &coeff) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4LivermorePhotoElectricModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 115 of file G4LivermorePhotoElectricModel.cc.

117 {
118  if (verboseLevel > 2) {
119  G4cout << "Calling G4LivermorePhotoElectricModel::Initialise()" << G4endl;
120  }
121 
122  if(IsMaster()) {
123 
124  if(!fWater) {
125  fWater = G4Material::GetMaterial("G4_WATER", false);
126  if(fWater) { fWaterEnergyLimit = 13.6*eV; }
127  }
128 
129  if(!fShellCrossSection) { fShellCrossSection = new G4ElementData(); }
130 
131  char* path = getenv("G4LEDATA");
132 
133  G4ProductionCutsTable* theCoupleTable =
135  G4int numOfCouples = theCoupleTable->GetTableSize();
136 
137  for(G4int i=0; i<numOfCouples; ++i) {
138  const G4MaterialCutsCouple* couple =
139  theCoupleTable->GetMaterialCutsCouple(i);
140  const G4Material* material = couple->GetMaterial();
141  const G4ElementVector* theElementVector = material->GetElementVector();
142  G4int nelm = material->GetNumberOfElements();
143 
144  for (G4int j=0; j<nelm; ++j) {
145  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
146  if(Z < 1) { Z = 1; }
147  else if(Z > maxZ) { Z = maxZ; }
148  if(!fCrossSection[Z]) { ReadData(Z, path); }
149  }
150  }
151  }
152  //
153  if (verboseLevel > 2) {
154  G4cout << "Loaded cross section files for LivermorePhotoElectric model"
155  << G4endl;
156  }
157  if(!isInitialised) {
158  isInitialised = true;
160 
161  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
162  }
163  fDeexcitationActive = false;
164  if(fAtomDeexcitation) {
165  fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
166  }
167 
168  if (verboseLevel > 0) {
169  G4cout << "LivermorePhotoElectric model is initialized " << G4endl
170  << G4endl;
171  }
172 }
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::vector< G4Element * > G4ElementVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * fParticleChange
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4LivermorePhotoElectricModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 596 of file G4LivermorePhotoElectricModel.cc.

598 {
599  G4AutoLock l(&LivermorePhotoElectricModelMutex);
600  // G4cout << "G4LivermorePhotoElectricModel::InitialiseForElement Z= "
601  // << Z << G4endl;
602  if(!fCrossSection[Z]) { ReadData(Z); }
603  l.unlock();
604 }

Here is the call graph for this function:

Here is the caller graph for this function:

void G4LivermorePhotoElectricModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 257 of file G4LivermorePhotoElectricModel.cc.

263 {
264  G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
265  if (verboseLevel > 3) {
266  G4cout << "G4LivermorePhotoElectricModel::SampleSecondaries() Egamma(keV)= "
267  << gammaEnergy/keV << G4endl;
268  }
269 
270  // kill incident photon
273 
274  // low-energy photo-effect in water - full absorption
275  const G4Material* material = couple->GetMaterial();
276  if(fWater && (material == fWater ||
277  material->GetBaseMaterial() == fWater)) {
278  if(gammaEnergy <= fWaterEnergyLimit) {
280  return;
281  }
282  }
283 
284  // Returns the normalized direction of the momentum
285  G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
286 
287  // Select randomly one element in the current material
288  //G4cout << "Select random atom Egamma(keV)= " << gammaEnergy/keV << G4endl;
289  const G4Element* elm = SelectRandomAtom(material, theGamma, gammaEnergy);
290  G4int Z = G4lrint(elm->GetZ());
291 
292  // Select the ionised shell in the current atom according to shell
293  // cross sections
294  // G4cout << "Select random shell Z= " << Z << G4endl;
295 
296  if(Z >= maxZ) { Z = maxZ-1; }
297 
298  // element was not initialised gamma should be absorbed
299  if(!fCrossSection[Z]) {
301  return;
302  }
303 
304  // shell index
305  size_t shellIdx = 0;
306  size_t nn = fNShellsUsed[Z];
307 
308  if(nn > 1) {
309  if(gammaEnergy >= (*(fParam[Z]))[0]) {
310  G4double x1 = 1.0/gammaEnergy;
311  G4double x2 = x1*x1;
312  G4double x3 = x2*x1;
313  G4double x4 = x3*x1;
314  G4int idx = nn*6 - 4;
315  // when do sampling common factors are not taken into account
316  // so cross section is not real
317  G4double cs0 = G4UniformRand()*((*(fParam[Z]))[idx]
318  + x1*(*(fParam[Z]))[idx+1]
319  + x2*(*(fParam[Z]))[idx+2]
320  + x3*(*(fParam[Z]))[idx+3]
321  + x4*(*(fParam[Z]))[idx+4]);
322  for(shellIdx=0; shellIdx<nn; ++shellIdx) {
323  idx = shellIdx*6 + 2;
324  if(gammaEnergy > (*(fParam[Z]))[idx-1]) {
325  G4double cs = (*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
326  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
327  + x4*(*(fParam[Z]))[idx+4];
328  if(cs >= cs0) { break; }
329  }
330  }
331  if(shellIdx >= nn) { shellIdx = nn-1; }
332 
333  } else {
334 
335  // when do sampling common factors are not taken into account
336  // so cross section is not real
337  G4double cs = G4UniformRand();
338 
339  if(gammaEnergy >= (*(fParam[Z]))[1]) {
340  cs *= (fCrossSection[Z])->Value(gammaEnergy);
341  } else {
342  cs *= (fCrossSectionLE[Z])->Value(gammaEnergy);
343  }
344 
345  for(size_t j=0; j<nn; ++j) {
346  shellIdx = (size_t)fShellCrossSection->GetComponentID(Z, j);
347  if(gammaEnergy > (*(fParam[Z]))[6*shellIdx+1]) {
348  cs -= fShellCrossSection->GetValueForComponent(Z, j, gammaEnergy);
349  }
350  if(cs <= 0.0 || j+1 == nn) { break; }
351  }
352  }
353  }
354 
355  G4double bindingEnergy = (*(fParam[Z]))[shellIdx*6 + 1];
356  //G4cout << "Z= " << Z << " shellIdx= " << shellIdx
357  // << " nShells= " << fNShells[Z]
358  // << " Ebind(keV)= " << bindingEnergy/keV
359  // << " Egamma(keV)= " << gammaEnergy/keV << G4endl;
360 
361  const G4AtomicShell* shell = 0;
362 
363  // no de-excitation from the last shell
364  if(fDeexcitationActive && shellIdx + 1 < nn) {
366  shell = fAtomDeexcitation->GetAtomicShell(Z, as);
367  }
368 
369  // If binding energy of the selected shell is larger than photon energy
370  // do not generate secondaries
371  if(gammaEnergy < bindingEnergy) {
373  return;
374  }
375 
376  // Primary outcoming electron
377  G4double eKineticEnergy = gammaEnergy - bindingEnergy;
378  G4double edep = bindingEnergy;
379 
380  // Calculate direction of the photoelectron
381  G4ThreeVector electronDirection =
382  GetAngularDistribution()->SampleDirection(aDynamicGamma,
383  eKineticEnergy,
384  shellIdx,
385  couple->GetMaterial());
386 
387  // The electron is created
388  G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
389  electronDirection,
390  eKineticEnergy);
391  fvect->push_back(electron);
392 
393  // Sample deexcitation
394  if(shell) {
395  G4int index = couple->GetIndex();
396  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
397  G4int nbefore = fvect->size();
398 
399  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
400  G4int nafter = fvect->size();
401  if(nafter > nbefore) {
402  G4double esec = 0.0;
403  for (G4int j=nbefore; j<nafter; ++j) {
404 
405  G4double e = ((*fvect)[j])->GetKineticEnergy();
406  if(esec + e > edep) {
407  // correct energy in order to have energy balance
408  e = edep - esec;
409  ((*fvect)[j])->SetKineticEnergy(e);
410  esec += e;
411  // delete the rest of secondaries (should not happens)
412  for (G4int jj=nafter-1; jj>j; --jj) {
413  delete (*fvect)[jj];
414  fvect->pop_back();
415  }
416  break;
417  }
418  esec += e;
419  }
420  edep -= esec;
421  }
422  }
423  }
424  // energy balance - excitation energy left
425  if(edep > 0.0) {
427  }
428 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
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:97
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)
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
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)
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * fParticleChange
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4LivermorePhotoElectricModel::SetLimitNumberOfShells ( G4int  n)
inline

Definition at line 122 of file G4LivermorePhotoElectricModel.hh.

123 {
124  nShellLimit = n;
125 }

Member Data Documentation

G4ParticleChangeForGamma* G4LivermorePhotoElectricModel::fParticleChange
protected

Definition at line 88 of file G4LivermorePhotoElectricModel.hh.


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