Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermorePolarizedPhotoElectricGDModel Class Reference

#include <G4LivermorePolarizedPhotoElectricGDModel.hh>

+ Inheritance 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

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

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 42 of file G4LivermorePolarizedPhotoElectricGDModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedPhotoElectricGDModel()

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

Definition at line 59 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

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

◆ ~G4LivermorePolarizedPhotoElectricGDModel()

G4LivermorePolarizedPhotoElectricGDModel::~G4LivermorePolarizedPhotoElectricGDModel ( )
virtual

Definition at line 89 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

90{
91 if(IsMaster()) {
92 delete fShellCrossSection;
93 for(G4int i=0; i<maxZ; ++i) {
94 delete fParam[i];
95 fParam[i] = 0;
96 delete fCrossSection[i];
97 fCrossSection[i] = 0;
98 delete fCrossSectionLE[i];
99 fCrossSectionLE[i] = 0;
100 }
101 }
102}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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 169 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 107 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

110{
111 if (verboseLevel > 2) {
112 G4cout << "Calling G4LivermorePolarizedPhotoElectricGDModel::Initialise()" << G4endl;
113 }
114
115 if(IsMaster()) {
116
117 if(!fWater) {
118 fWater = G4Material::GetMaterial("G4_WATER", false);
119 if(fWater) { fWaterEnergyLimit = 13.6*eV; }
120 }
121
122 if(!fShellCrossSection) { fShellCrossSection = new G4ElementData(); }
123
124 char* path = std::getenv("G4LEDATA");
125
126 G4ProductionCutsTable* theCoupleTable =
128 G4int numOfCouples = theCoupleTable->GetTableSize();
129
130 for(G4int i=0; i<numOfCouples; ++i) {
131 const G4MaterialCutsCouple* couple =
132 theCoupleTable->GetMaterialCutsCouple(i);
133 const G4Material* material = couple->GetMaterial();
134 const G4ElementVector* theElementVector = material->GetElementVector();
135 G4int nelm = material->GetNumberOfElements();
136
137 for (G4int j=0; j<nelm; ++j) {
138 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
139 if(Z < 1) { Z = 1; }
140 else if(Z > maxZ) { Z = maxZ; }
141 if(!fCrossSection[Z]) { ReadData(Z, path); }
142 }
143 }
144 }
145 //
146 if (verboseLevel > 2) {
147 G4cout << "Loaded cross section files for LivermorePhotoElectric model"
148 << G4endl;
149 }
150 if(!isInitialised) {
151 isInitialised = true;
153
154 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
155 }
156 fDeexcitationActive = false;
157 if(fAtomDeexcitation) {
158 fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
159 }
160
161 if (verboseLevel > 0) {
162 G4cout << "LivermorePolarizedPhotoElectric model is initialized " << G4endl
163 << G4endl;
164 }
165}
std::vector< G4Element * > G4ElementVector
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 744 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

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

Referenced by ComputeCrossSectionPerAtom().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 221 of file G4LivermorePolarizedPhotoElectricGDModel.cc.

227{
228 if (verboseLevel > 3) {
229 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedPhotoElectricGDModel"
230 << G4endl;
231 }
232
233 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
234 if (verboseLevel > 3) {
235 G4cout << "G4LivermorePolarizedPhotoElectricGDModel::SampleSecondaries() Egamma(keV)= "
236 << photonEnergy/keV << G4endl;
237 }
238
239 G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
240 G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
241
242 // kill incident photon
245
246 // low-energy photo-effect in water - full absorption
247
248 const G4Material* material = couple->GetMaterial();
249 if(fWater && (material == fWater ||
250 material->GetBaseMaterial() == fWater)) {
251 if(photonEnergy <= fWaterEnergyLimit) {
253 return;
254 }
255 }
256
257 // Protection: a polarisation parallel to the
258 // direction causes problems;
259 // in that case find a random polarization
260
261 // Make sure that the polarization vector is perpendicular to the
262 // gamma direction. If not
263
264 if(!(gammaPolarization0.isOrthogonal(photonDirection, 1e-6))||(gammaPolarization0.mag()==0))
265 { // only for testing now
266 gammaPolarization0 = GetRandomPolarization(photonDirection);
267 }
268 else
269 {
270 if ( gammaPolarization0.howOrthogonal(photonDirection) != 0)
271 {
272 gammaPolarization0 = GetPerpendicularPolarization(photonDirection, gammaPolarization0);
273 }
274 }
275
276 // End of Protection
277
278 // G4double E0_m = photonEnergy / electron_mass_c2 ;
279
280 // Shell
281
282 // Select randomly one element in the current material
283 //G4cout << "Select random atom Egamma(keV)= " << photonEnergy/keV << G4endl;
284 const G4Element* elm = SelectRandomAtom(material, theGamma, photonEnergy);
285 G4int Z = G4lrint(elm->GetZ());
286
287
288 // Select the ionised shell in the current atom according to shell
289 // cross sections
290 // G4cout << "Select random shell Z= " << Z << G4endl;
291
292 if(Z >= maxZ) { Z = maxZ-1; }
293
294 // element was not initialised gamma should be absorbed
295 if(!fCrossSection[Z]) {
297 return;
298 }
299
300 // shell index
301 size_t shellIdx = 0;
302 size_t nn = fNShellsUsed[Z];
303
304 if(nn > 1) {
305 if(photonEnergy >= (*(fParam[Z]))[0]) {
306 G4double x1 = 1.0/photonEnergy;
307 G4double x2 = x1*x1;
308 G4double x3 = x2*x1;
309 G4double x4 = x3*x1;
310 G4int idx = nn*6 - 4;
311 // when do sampling common factors are not taken into account
312 // so cross section is not real
313 G4double cs0 = G4UniformRand()*((*(fParam[Z]))[idx]
314 + x1*(*(fParam[Z]))[idx+1]
315 + x2*(*(fParam[Z]))[idx+2]
316 + x3*(*(fParam[Z]))[idx+3]
317 + x4*(*(fParam[Z]))[idx+4]);
318 for(shellIdx=0; shellIdx<nn; ++shellIdx) {
319 idx = shellIdx*6 + 2;
320 if(photonEnergy > (*(fParam[Z]))[idx-1]) {
321 G4double cs = (*(fParam[Z]))[idx] + x1*(*(fParam[Z]))[idx+1]
322 + x2*(*(fParam[Z]))[idx+2] + x3*(*(fParam[Z]))[idx+3]
323 + x4*(*(fParam[Z]))[idx+4];
324 if(cs >= cs0) { break; }
325 }
326 }
327 if(shellIdx >= nn) { shellIdx = nn-1; }
328
329 } else {
330
331 // when do sampling common factors are not taken into account
332 // so cross section is not real
333 G4double cs = G4UniformRand();
334
335 if(photonEnergy >= (*(fParam[Z]))[1]) {
336 cs *= (fCrossSection[Z])->Value(photonEnergy);
337 } else {
338 cs *= (fCrossSectionLE[Z])->Value(photonEnergy);
339 }
340
341 for(size_t j=0; j<nn; ++j) {
342 shellIdx = (size_t)fShellCrossSection->GetComponentID(Z, j);
343 if(photonEnergy > (*(fParam[Z]))[6*shellIdx+1]) {
344 cs -= fShellCrossSection->GetValueForComponent(Z, j, photonEnergy);
345 }
346 if(cs <= 0.0 || j+1 == nn) { break; }
347 }
348 }
349 }
350
351 G4double bindingEnergy = (*(fParam[Z]))[shellIdx*6 + 1];
352 //G4cout << "Z= " << Z << " shellIdx= " << shellIdx
353 // << " nShells= " << fNShells[Z]
354 // << " Ebind(keV)= " << bindingEnergy/keV
355 // << " Egamma(keV)= " << photonEnergy/keV << G4endl;
356
357 const G4AtomicShell* shell = 0;
358
359 // no de-excitation from the last shell
360 if(fDeexcitationActive && shellIdx + 1 < nn) {
362 shell = fAtomDeexcitation->GetAtomicShell(Z, as);
363 }
364
365 // If binding energy of the selected shell is larger than photon energy
366 // do not generate secondaries
367 if(photonEnergy < bindingEnergy) {
369 return;
370 }
371
372
373 // Electron
374
375 G4double eKineticEnergy = photonEnergy - bindingEnergy;
376 G4double edep = bindingEnergy;
377
378 G4double costheta = SetCosTheta(eKineticEnergy);
379 G4double sintheta = sqrt(1. - costheta*costheta);
380 G4double phi = SetPhi(photonEnergy,eKineticEnergy,costheta);
381 G4double dirX = sintheta*cos(phi);
382 G4double dirY = sintheta*sin(phi);
383 G4double dirZ = costheta;
384 G4ThreeVector electronDirection(dirX, dirY, dirZ);
385 SystemOfRefChange(photonDirection, electronDirection, gammaPolarization0);
387 electronDirection,
388 eKineticEnergy);
389 fvect->push_back(electron);
390
391 // Deexcitation
392 // Sample deexcitation
393 if(shell) {
394 G4int index = couple->GetIndex();
395 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
396 G4int nbefore = fvect->size();
397
398 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
399 G4int nafter = fvect->size();
400 if(nafter > nbefore) {
401 G4double esec = 0.0;
402 for (G4int j=nbefore; j<nafter; ++j) {
403
404 G4double e = ((*fvect)[j])->GetKineticEnergy();
405 if(esec + e > edep) {
406 // correct energy in order to have energy balance
407 e = edep - esec;
408 ((*fvect)[j])->SetKineticEnergy(e);
409 esec += e;
410 // delete the rest of secondaries (should not happens)
411 for (G4int jj=nafter-1; jj>j; --jj) {
412 delete (*fvect)[jj];
413 fvect->pop_back();
414 }
415 break;
416 }
417 esec += e;
418 }
419 edep -= esec;
420 }
421 }
422 }
423 // energy balance - excitation energy left
424 if(edep > 0.0) {
426 }
427}
G4AtomicShellEnumerator
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:233
double mag() const
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:215
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4int GetComponentID(G4int Z, size_t idx)
G4double GetValueForComponent(G4int Z, size_t idx, G4double kinEnergy)
G4double GetZ() const
Definition: G4Element.hh:130
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:231
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SetLimitNumberOfShells()

void G4LivermorePolarizedPhotoElectricGDModel::SetLimitNumberOfShells ( G4int  n)
inline

Definition at line 123 of file G4LivermorePolarizedPhotoElectricGDModel.hh.

124{
125 nShellLimit = n;
126}

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricGDModel::fParticleChange
protected

Definition at line 74 of file G4LivermorePolarizedPhotoElectricGDModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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