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

#include <G4LivermorePolarizedPhotoElectricGDModel.hh>

Inheritance diagram for G4LivermorePolarizedPhotoElectricGDModel:
Collaboration diagram for G4LivermorePolarizedPhotoElectricGDModel:

Public Member Functions

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

Constructor & Destructor Documentation

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

Definition at line 60 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

62  :G4VEmModel(nam),fParticleChange(nullptr),maxZ(99),
63  nShellLimit(100), fDeexcitationActive(false), isInitialised(false),
64  fAtomDeexcitation(nullptr)
65 {
66  verboseLevel= 0;
67  // Verbosity scale:
68  // 0 = nothing
69  // 1 = warning for energy non-conservation
70  // 2 = details of energy budget
71  // 3 = calculation of cross sections, file openings, sampling of atoms
72  // 4 = entering in methods
73 
74  theGamma = G4Gamma::Gamma();
75  theElectron = G4Electron::Electron();
76 
77  SetDeexcitationFlag(true);
78  fSandiaCof.resize(4,0.0);
79  fCurrSection = 0.0;
80 
81  if (verboseLevel > 0) {
82  G4cout << "Livermore Polarized PhotoElectric is constructed "
83  << " nShellLimit "
84  << nShellLimit << G4endl;
85  }
86 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
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:

G4LivermorePolarizedPhotoElectricGDModel::~G4LivermorePolarizedPhotoElectricGDModel ( )
virtual

Definition at line 90 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

91 {
92  if(IsMaster()) {
93  delete fShellCrossSection;
94  for(G4int i=0; i<maxZ; ++i) {
95  delete fParam[i];
96  fParam[i] = 0;
97  delete fCrossSection[i];
98  fCrossSection[i] = 0;
99  delete fCrossSectionLE[i];
100  fCrossSectionLE[i] = 0;
101  }
102  }
103 }
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 G4LivermorePolarizedPhotoElectricGDModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 170 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

175 {
176  if (verboseLevel > 3) {
177  G4cout << "G4LivermorePolarizedPhotoElectricGDModel::ComputeCrossSectionPerAtom():"
178  << " Z= " << ZZ << " R(keV)= " << GammaEnergy/keV << G4endl;
179  }
180  G4double cs = 0.0;
181  G4int Z = G4lrint(ZZ);
182  if(Z < 1 || Z >= maxZ) { return cs; }
183 
184  // if element was not initialised
185  // do initialisation safely for MT mode
186  if(!fCrossSection[Z]) {
187  InitialiseForElement(0, Z);
188  if(!fCrossSection[Z]) { return cs; }
189  }
190 
191  G4int idx = fNShells[Z]*6 - 4;
192  if (GammaEnergy < (*(fParam[Z]))[idx-1]) { GammaEnergy = (*(fParam[Z]))[idx-1]; }
193 
194  G4double x1 = 1.0/GammaEnergy;
195  G4double x2 = x1*x1;
196  G4double x3 = x2*x1;
197 
198  // parameterisation
199  if(GammaEnergy >= (*(fParam[Z]))[0]) {
200  G4double x4 = x2*x2;
201  cs = x1*((*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
202  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
203  + x4*(*(fParam[Z]))[idx+4]);
204  // high energy part
205  } else if (GammaEnergy >= (*(fParam[Z]))[1]) {
206  cs = x3*(fCrossSection[Z])->Value(GammaEnergy);
207 
208  // low energy part
209  } else {
210  cs = x3*(fCrossSectionLE[Z])->Value(GammaEnergy);
211  }
212  if (verboseLevel > 1) {
213  G4cout << "LivermorePolarizedPhotoElectricGDModel: E(keV)= " << GammaEnergy/keV
214  << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
215  }
216  return cs;
217 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
int G4lrint(double ad)
Definition: templates.hh:163
#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:

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

Implements G4VEmModel.

Definition at line 108 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

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

Reimplemented from G4VEmModel.

Definition at line 745 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

747 {
748  G4AutoLock l(&LivermorePolarizedPhotoElectricGDModelMutex);
749  // G4cout << "G4LivermorePhotoElectricModel::InitialiseForElement Z= "
750  // << Z << G4endl;
751  if(!fCrossSection[Z]) { ReadData(Z); }
752  l.unlock();
753 }

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 222 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

228 {
229  if (verboseLevel > 3) {
230  G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedPhotoElectricGDModel"
231  << G4endl;
232  }
233 
234  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
235  if (verboseLevel > 3) {
236  G4cout << "G4LivermorePolarizedPhotoElectricGDModel::SampleSecondaries() Egamma(keV)= "
237  << photonEnergy/keV << G4endl;
238  }
239 
240  G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
241  G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
242 
243  // kill incident photon
246 
247  // low-energy photo-effect in water - full absorption
248 
249  const G4Material* material = couple->GetMaterial();
250  if(fWater && (material == fWater ||
251  material->GetBaseMaterial() == fWater)) {
252  if(photonEnergy <= fWaterEnergyLimit) {
254  return;
255  }
256  }
257 
258  // Protection: a polarisation parallel to the
259  // direction causes problems;
260  // in that case find a random polarization
261 
262  // Make sure that the polarization vector is perpendicular to the
263  // gamma direction. If not
264 
265  if(!(gammaPolarization0.isOrthogonal(photonDirection, 1e-6))||(gammaPolarization0.mag()==0))
266  { // only for testing now
267  gammaPolarization0 = GetRandomPolarization(photonDirection);
268  }
269  else
270  {
271  if ( gammaPolarization0.howOrthogonal(photonDirection) != 0)
272  {
273  gammaPolarization0 = GetPerpendicularPolarization(photonDirection, gammaPolarization0);
274  }
275  }
276 
277  // End of Protection
278 
279  // G4double E0_m = photonEnergy / electron_mass_c2 ;
280 
281  // Shell
282 
283  // Select randomly one element in the current material
284  //G4cout << "Select random atom Egamma(keV)= " << photonEnergy/keV << G4endl;
285  const G4Element* elm = SelectRandomAtom(material, theGamma, photonEnergy);
286  G4int Z = G4lrint(elm->GetZ());
287 
288 
289  // Select the ionised shell in the current atom according to shell
290  // cross sections
291  // G4cout << "Select random shell Z= " << Z << G4endl;
292 
293  if(Z >= maxZ) { Z = maxZ-1; }
294 
295  // element was not initialised gamma should be absorbed
296  if(!fCrossSection[Z]) {
298  return;
299  }
300 
301  // shell index
302  size_t shellIdx = 0;
303  size_t nn = fNShellsUsed[Z];
304 
305  if(nn > 1) {
306  if(photonEnergy >= (*(fParam[Z]))[0]) {
307  G4double x1 = 1.0/photonEnergy;
308  G4double x2 = x1*x1;
309  G4double x3 = x2*x1;
310  G4double x4 = x3*x1;
311  G4int idx = nn*6 - 4;
312  // when do sampling common factors are not taken into account
313  // so cross section is not real
314  G4double cs0 = G4UniformRand()*((*(fParam[Z]))[idx]
315  + x1*(*(fParam[Z]))[idx+1]
316  + x2*(*(fParam[Z]))[idx+2]
317  + x3*(*(fParam[Z]))[idx+3]
318  + x4*(*(fParam[Z]))[idx+4]);
319  for(shellIdx=0; shellIdx<nn; ++shellIdx) {
320  idx = shellIdx*6 + 2;
321  if(photonEnergy > (*(fParam[Z]))[idx-1]) {
322  G4double cs = (*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
323  + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
324  + x4*(*(fParam[Z]))[idx+4];
325  if(cs >= cs0) { break; }
326  }
327  }
328  if(shellIdx >= nn) { shellIdx = nn-1; }
329 
330  } else {
331 
332  // when do sampling common factors are not taken into account
333  // so cross section is not real
334  G4double cs = G4UniformRand();
335 
336  if(photonEnergy >= (*(fParam[Z]))[1]) {
337  cs *= (fCrossSection[Z])->Value(photonEnergy);
338  } else {
339  cs *= (fCrossSectionLE[Z])->Value(photonEnergy);
340  }
341 
342  for(size_t j=0; j<nn; ++j) {
343  shellIdx = (size_t)fShellCrossSection->GetComponentID(Z, j);
344  if(photonEnergy > (*(fParam[Z]))[6*shellIdx+1]) {
345  cs -= fShellCrossSection->GetValueForComponent(Z, j, photonEnergy);
346  }
347  if(cs <= 0.0 || j+1 == nn) { break; }
348  }
349  }
350  }
351 
352  G4double bindingEnergy = (*(fParam[Z]))[shellIdx*6 + 1];
353  //G4cout << "Z= " << Z << " shellIdx= " << shellIdx
354  // << " nShells= " << fNShells[Z]
355  // << " Ebind(keV)= " << bindingEnergy/keV
356  // << " Egamma(keV)= " << photonEnergy/keV << G4endl;
357 
358  const G4AtomicShell* shell = 0;
359 
360  // no de-excitation from the last shell
361  if(fDeexcitationActive && shellIdx + 1 < nn) {
363  shell = fAtomDeexcitation->GetAtomicShell(Z, as);
364  }
365 
366  // If binding energy of the selected shell is larger than photon energy
367  // do not generate secondaries
368  if(photonEnergy < bindingEnergy) {
370  return;
371  }
372 
373 
374  // Electron
375 
376  G4double eKineticEnergy = photonEnergy - bindingEnergy;
377  G4double edep = bindingEnergy;
378 
379  G4double costheta = SetCosTheta(eKineticEnergy);
380  G4double sintheta = sqrt(1. - costheta*costheta);
381  G4double phi = SetPhi(photonEnergy,eKineticEnergy,costheta);
382  G4double dirX = sintheta*cos(phi);
383  G4double dirY = sintheta*sin(phi);
384  G4double dirZ = costheta;
385  G4ThreeVector electronDirection(dirX, dirY, dirZ);
386  SystemOfRefChange(photonDirection, electronDirection, gammaPolarization0);
388  electronDirection,
389  eKineticEnergy);
390  fvect->push_back(electron);
391 
392  // Deexcitation
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
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
int G4int
Definition: G4Types.hh:78
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
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
const G4ThreeVector & GetMomentumDirection() const
int G4lrint(double ad)
Definition: templates.hh:163
const G4ThreeVector & GetPolarization() const
G4double GetValueForComponent(G4int Z, size_t idx, G4double kinEnergy)
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
static G4Electron * Electron()
Definition: G4Electron.cc:94
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)
double mag() const
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 G4LivermorePolarizedPhotoElectricGDModel::SetLimitNumberOfShells ( G4int  n)
inline

Definition at line 124 of file G4LivermorePolarizedPhotoElectricGDModel.hh.

125 {
126  nShellLimit = n;
127 }

Member Data Documentation

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricGDModel::fParticleChange
protected

Definition at line 75 of file G4LivermorePolarizedPhotoElectricGDModel.hh.


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