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

#include <G4LivermorePolarizedPhotoElectricModel.hh>

+ Inheritance 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 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 49 of file G4LivermorePolarizedPhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedPhotoElectricModel()

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

Definition at line 64 of file G4LivermorePolarizedPhotoElectricModel.cc.

65 : G4VEmModel(nam),fParticleChange(nullptr),maxZ(99),
66 nShellLimit(100),fDeexcitationActive(false),isInitialised(false),
67 fAtomDeexcitation(nullptr)
68{
69 verboseLevel= 0;
70 // Verbosity scale:
71 // 0 = nothing
72 // 1 = warning for energy non-conservation
73 // 2 = details of energy budget
74 // 3 = calculation of cross sections, file openings, sampling of atoms
75 // 4 = entering in methods
76
77 theGamma = G4Gamma::Gamma();
78 theElectron = G4Electron::Electron();
79
80 // default generator
81 // SetAngularDistribution(new G4SauterGavrilaAngularDistribution());
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
91 fSandiaCof.resize(4,0.0);
92 fCurrSection = 0.0;
93
94}
#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
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4LivermorePolarizedPhotoElectricModel()

G4LivermorePolarizedPhotoElectricModel::~G4LivermorePolarizedPhotoElectricModel ( )
virtual

Definition at line 98 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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 206 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

◆ CrossSectionPerVolume()

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 177 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 116 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

Reimplemented from G4VEmModel.

Definition at line 597 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

Referenced by ComputeCrossSectionPerAtom().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 258 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

Definition at line 121 of file G4LivermorePolarizedPhotoElectricModel.hh.

122{
123 nShellLimit = n;
124}

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricModel::fParticleChange
protected

Definition at line 87 of file G4LivermorePolarizedPhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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