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

#include <G4LivermorePolarizedPhotoElectricModel.hh>

Inheritance diagram for G4LivermorePolarizedPhotoElectricModel:
Collaboration diagram for G4LivermorePolarizedPhotoElectricModel:

Public Member Functions

 G4LivermorePolarizedPhotoElectricModel (const G4String &nam="LivermorePolarizedPhotoElectric")
 
virtual ~G4LivermorePolarizedPhotoElectricModel ()
 
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 G4LivermorePolarizedPhotoElectricModel.hh.

Constructor & Destructor Documentation

G4LivermorePolarizedPhotoElectricModel::G4LivermorePolarizedPhotoElectricModel ( const G4String nam = "LivermorePolarizedPhotoElectric")

Definition at line 65 of file G4LivermorePolarizedPhotoElectricModel.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
82  // SetAngularDistribution(new G4SauterGavrilaAngularDistribution());
84 
85  if(verboseLevel>0) {
86  G4cout << "Livermore PhotoElectric is constructed "
87  << " nShellLimit= " << nShellLimit << G4endl;
88  }
89 
90  //Mark this model as "applicable" for atomic deexcitation
91  SetDeexcitationFlag(true);
92  fSandiaCof.resize(4,0.0);
93  fCurrSection = 0.0;
94 
95 }
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

Here is the call graph for this function:

G4LivermorePolarizedPhotoElectricModel::~G4LivermorePolarizedPhotoElectricModel ( )
virtual

Definition at line 99 of file G4LivermorePolarizedPhotoElectricModel.cc.

100 {
101  if(IsMaster()) {
102  delete fShellCrossSection;
103  for(G4int i=0; i<maxZ; ++i) {
104  delete fParam[i];
105  fParam[i] = 0;
106  delete fCrossSection[i];
107  fCrossSection[i] = 0;
108  delete fCrossSectionLE[i];
109  fCrossSectionLE[i] = 0;
110  }
111  }
112 }
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 G4LivermorePolarizedPhotoElectricModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  energy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 207 of file G4LivermorePolarizedPhotoElectricModel.cc.

212 {
213  if (verboseLevel > 3) {
214  G4cout << "G4LivermorePolarizedPhotoElectricModel::ComputeCrossSectionPerAtom():"
215  << " Z= " << ZZ << " R(keV)= " << energy/keV << G4endl;
216  }
217  G4double cs = 0.0;
218  G4int Z = G4lrint(ZZ);
219  if(Z < 1 || Z >= maxZ) { return cs; }
220 
221  // if element was not initialised
222  // do initialisation safely for MT mode
223  if(!fCrossSection[Z]) {
224  InitialiseForElement(0, Z);
225  if(!fCrossSection[Z]) { return cs; }
226  }
227 
228  G4int idx = fNShells[Z]*6 - 4;
229  if (energy < (*(fParam[Z]))[idx-1]) { energy = (*(fParam[Z]))[idx-1]; }
230 
231  G4double x1 = 1.0/energy;
232  G4double x2 = x1*x1;
233  G4double x3 = x2*x1;
234 
235  // parameterisation
236  if(energy >= (*(fParam[Z]))[0]) {
237  G4double x4 = x2*x2;
238  cs = x1*((*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
239  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
240  + x4*(*(fParam[Z]))[idx+4]);
241  // high energy part
242  } else if(energy >= (*(fParam[Z]))[1]) {
243  cs = x3*(fCrossSection[Z])->Value(energy);
244 
245  // low energy part
246  } else {
247  cs = x3*(fCrossSectionLE[Z])->Value(energy);
248  }
249  if (verboseLevel > 1) {
250  G4cout << "LivermorePolarizedPhotoElectricModel: E(keV)= " << energy/keV
251  << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
252  }
253  return cs;
254 }
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 G4LivermorePolarizedPhotoElectricModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  energy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 178 of file G4LivermorePolarizedPhotoElectricModel.cc.

183 {
184  fCurrSection = 0.0;
185  if(fWater && (material == fWater ||
186  material->GetBaseMaterial() == fWater)) {
187  if(energy <= fWaterEnergyLimit) {
188  fWater->GetSandiaTable()->GetSandiaCofWater(energy, fSandiaCof);
189 
190  G4double energy2 = energy*energy;
191  G4double energy3 = energy*energy2;
192  G4double energy4 = energy2*energy2;
193 
194  fCurrSection = material->GetDensity()*
195  (fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
196  fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4);
197  }
198  }
199  if(0.0 == fCurrSection) {
200  fCurrSection = G4VEmModel::CrossSectionPerVolume(material, p, energy);
201  }
202  return fCurrSection;
203 }
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 G4LivermorePolarizedPhotoElectricModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 117 of file G4LivermorePolarizedPhotoElectricModel.cc.

119 {
120  if (verboseLevel > 2) {
121  G4cout << "Calling G4LivermorePolarizedPhotoElectricModel::Initialise()" << G4endl;
122  }
123 
124  if(IsMaster()) {
125 
126  if(!fWater) {
127  fWater = G4Material::GetMaterial("G4_WATER", false);
128  if(fWater) { fWaterEnergyLimit = 13.6*eV; }
129  }
130 
131  if(!fShellCrossSection) { fShellCrossSection = new G4ElementData(); }
132 
133  char* path = getenv("G4LEDATA");
134 
135  G4ProductionCutsTable* theCoupleTable =
137  G4int numOfCouples = theCoupleTable->GetTableSize();
138 
139  for(G4int i=0; i<numOfCouples; ++i) {
140  const G4MaterialCutsCouple* couple =
141  theCoupleTable->GetMaterialCutsCouple(i);
142  const G4Material* material = couple->GetMaterial();
143  const G4ElementVector* theElementVector = material->GetElementVector();
144  G4int nelm = material->GetNumberOfElements();
145 
146  for (G4int j=0; j<nelm; ++j) {
147  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
148  if(Z < 1) { Z = 1; }
149  else if(Z > maxZ) { Z = maxZ; }
150  if(!fCrossSection[Z]) { ReadData(Z, path); }
151  }
152  }
153  }
154  //
155  if (verboseLevel > 2) {
156  G4cout << "Loaded cross section files for LivermorePolarizedPhotoElectric model"
157  << G4endl;
158  }
159  if(!isInitialised) {
160  isInitialised = true;
162 
163  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
164  }
165  fDeexcitationActive = false;
166  if(fAtomDeexcitation) {
167  fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
168  }
169 
170  if (verboseLevel > 0) {
171  G4cout << "LivermorePolarizedPhotoElectric model is initialized " << G4endl
172  << G4endl;
173  }
174 }
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 * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 598 of file G4LivermorePolarizedPhotoElectricModel.cc.

600 {
601  G4AutoLock l(&LivermorePolarizedPhotoElectricModelMutex);
602  // G4cout << "G4LivermorePolarizedPhotoElectricModel::InitialiseForElement Z= "
603  // << Z << G4endl;
604  if(!fCrossSection[Z]) { ReadData(Z); }
605  l.unlock();
606 }

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 259 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

Definition at line 122 of file G4LivermorePolarizedPhotoElectricModel.hh.

123 {
124  nShellLimit = n;
125 }

Member Data Documentation

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricModel::fParticleChange
protected

Definition at line 88 of file G4LivermorePolarizedPhotoElectricModel.hh.


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