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

#include <G4LivermorePhotoElectricModel.hh>

+ Inheritance diagram for G4LivermorePhotoElectricModel:

Public Member Functions

 G4LivermorePhotoElectricModel (const G4String &nam="LivermorePhElectric")
 
virtual ~G4LivermorePhotoElectricModel ()
 
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)
 
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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 50 of file G4LivermorePhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePhotoElectricModel()

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

Definition at line 54 of file G4LivermorePhotoElectricModel.cc.

56 : G4VEmModel(nam),fParticleChange(0),maxZ(99),
57 nShellLimit(100),fDeexcitationActive(false),isInitialised(false),
58 fAtomDeexcitation(0)
59{
60 verboseLevel= 0;
61 // Verbosity scale:
62 // 0 = nothing
63 // 1 = warning for energy non-conservation
64 // 2 = details of energy budget
65 // 3 = calculation of cross sections, file openings, sampling of atoms
66 // 4 = entering in methods
67
68 theGamma = G4Gamma::Gamma();
69 theElectron = G4Electron::Electron();
70
71 // default generator
73
74 for(G4int i=0; i<maxZ; ++i) {
75 fCrossSection[i] = 0;
76 fCrossSectionLE[i] = 0;
77 fNShells[i] = 0;
78 fNShellsUsed[i] = 0;
79 }
80
81 if(verboseLevel>0) {
82 G4cout << "Livermore PhotoElectric is constructed "
83 << " nShellLimit= " << nShellLimit << G4endl;
84 }
85
86 //Mark this model as "applicable" for atomic deexcitation
88}
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ParticleChangeForGamma * fParticleChange
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515

◆ ~G4LivermorePhotoElectricModel()

G4LivermorePhotoElectricModel::~G4LivermorePhotoElectricModel ( )
virtual

Definition at line 92 of file G4LivermorePhotoElectricModel.cc.

93{
94 for(G4int i=0; i<maxZ; ++i) {
95 delete fCrossSection[i];
96 delete fCrossSectionLE[i];
97 }
98}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 155 of file G4LivermorePhotoElectricModel.cc.

160{
161 if (verboseLevel > 3) {
162 G4cout << "G4LivermorePhotoElectricModel::Calling ComputeCrossSectionPerAtom()"
163 << " Z= " << ZZ << " R(keV)= " << energy/keV << G4endl;
164 }
165 G4double cs = 0.0;
166 G4double gammaEnergy = energy;
167
168 G4int Z = G4lrint(ZZ);
169 if(Z < 1 || Z >= maxZ) { return cs; }
170
171 // element was not initialised
172 if(!fCrossSection[Z]) {
173 char* path = getenv("G4LEDATA");
174 ReadData(Z, path);
175 if(!fCrossSection[Z]) { return cs; }
176 }
177
178 G4int idx = fNShells[Z]*6 - 4;
179 if (gammaEnergy <= (fParam[Z])[idx-1]) { return cs; }
180
181 G4double x1 = 1.0/gammaEnergy;
182 G4double x2 = x1*x1;
183 G4double x3 = x2*x1;
184
185 // parameterisation
186 if(gammaEnergy >= (fParam[Z])[0]) {
187 G4double x4 = x2*x2;
188 cs = x1*((fParam[Z])[idx] + x1*(fParam[Z])[idx+1]
189 + x2*(fParam[Z])[idx+2] + x3*(fParam[Z])[idx+3]
190 + x4*(fParam[Z])[idx+4]);
191 // high energy part
192 } else if(gammaEnergy >= (fParam[Z])[1]) {
193 cs = x3*(fCrossSection[Z])->Value(gammaEnergy);
194
195 // low energy part
196 } else {
197 cs = x3*(fCrossSectionLE[Z])->Value(gammaEnergy);
198 }
199 if (verboseLevel > 1) {
200 G4cout << "LivermorePhotoElectricModel: E(keV)= " << gammaEnergy/keV
201 << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
202 }
203 return cs;
204}
double G4double
Definition: G4Types.hh:64
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:286
int G4lrint(double ad)
Definition: templates.hh:163

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 103 of file G4LivermorePhotoElectricModel.cc.

105{
106 if (verboseLevel > 2) {
107 G4cout << "Calling G4LivermorePhotoElectricModel::Initialise()" << G4endl;
108 }
109
110 char* path = getenv("G4LEDATA");
111
112 G4ProductionCutsTable* theCoupleTable =
114 G4int numOfCouples = theCoupleTable->GetTableSize();
115
116 for(G4int i=0; i<numOfCouples; ++i)
117 {
118 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
119 const G4Material* material = couple->GetMaterial();
120 const G4ElementVector* theElementVector = material->GetElementVector();
121 G4int nelm = material->GetNumberOfElements();
122
123 for (G4int j=0; j<nelm; ++j)
124 {
125 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
126 if(Z < 1) { Z = 1; }
127 else if(Z > maxZ) { Z = maxZ; }
128 if(!fCrossSection[Z]) { ReadData(Z, path); }
129 }
130 }
131 //
132 if (verboseLevel > 2) {
133 G4cout << "Loaded cross section files for LivermorePhotoElectric model"
134 << G4endl;
135 }
136 if(!isInitialised) {
137 isInitialised = true;
139
140 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
141 }
142 fDeexcitationActive = false;
143 if(fAtomDeexcitation) {
144 fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
145 }
146
147 if (verboseLevel > 0) {
148 G4cout << "LivermorePhotoElectric model is initialized " << G4endl
149 << G4endl;
150 }
151}
std::vector< G4Element * > G4ElementVector
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 209 of file G4LivermorePhotoElectricModel.cc.

215{
216 G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
217 if (verboseLevel > 3) {
218 G4cout << "G4LivermorePhotoElectricModel::SampleSecondaries() Egamma(keV)= "
219 << gammaEnergy/keV << G4endl;
220 }
221
222 // kill incident photon
225
226 // Returns the normalized direction of the momentum
227 G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
228
229 // Select randomly one element in the current material
230 //G4cout << "Select random atom Egamma(keV)= " << gammaEnergy/keV << G4endl;
231 const G4Element* elm = SelectRandomAtom(couple->GetMaterial(),theGamma,
232 gammaEnergy);
233 G4int Z = G4lrint(elm->GetZ());
234
235 // Select the ionised shell in the current atom according to shell
236 // cross sections
237 // G4cout << "Select random shell Z= " << Z << G4endl;
238
239 if(Z >= maxZ) { Z = maxZ-1; }
240
241 // element was not initialised
242 if(!fCrossSection[Z]) {
243 char* path = getenv("G4LEDATA");
244 ReadData(Z, path);
245 if(!fCrossSection[Z]) {
247 return;
248 }
249 }
250
251 // shell index
252 size_t shellIdx = 0;
253 size_t nn = fNShellsUsed[Z];
254
255 if(nn > 1) {
256 if(gammaEnergy >= (fParam[Z])[0]) {
257 G4double x1 = 1.0/gammaEnergy;
258 G4double x2 = x1*x1;
259 G4double x3 = x2*x1;
260 G4double x4 = x3*x1;
261 G4int idx = nn*6 - 4;
262 // when do sampling common factors are not taken into account
263 // so cross section is not real
264 G4double cs0 = G4UniformRand()*((fParam[Z])[idx] + x1*(fParam[Z])[idx+1]
265 + x2*(fParam[Z])[idx+2]
266 + x3*(fParam[Z])[idx+3]
267 + x4*(fParam[Z])[idx+4]);
268 for(shellIdx=0; shellIdx<nn; ++shellIdx) {
269 idx = shellIdx*6 + 2;
270 if(gammaEnergy > (fParam[Z])[idx-1]) {
271 G4double cs = (fParam[Z])[idx] + x1*(fParam[Z])[idx+1]
272 + x2*(fParam[Z])[idx+2] + x3*(fParam[Z])[idx+3]
273 + x4*(fParam[Z])[idx+4];
274 if(cs >= cs0) { break; }
275 }
276 }
277 if(shellIdx >= nn) { shellIdx = nn-1; }
278
279 } else {
280
281 // when do sampling common factors are not taken into account
282 // so cross section is not real
283 G4double cs = G4UniformRand();
284
285 if(gammaEnergy >= (fParam[Z])[1]) {
286 cs *= (fCrossSection[Z])->Value(gammaEnergy);
287 } else {
288 cs *= (fCrossSectionLE[Z])->Value(gammaEnergy);
289 }
290
291 for(size_t j=0; j<nn; ++j) {
292 shellIdx = (size_t)fShellCrossSection.GetComponentID(Z, j);
293 if(gammaEnergy > (fParam[Z])[6*shellIdx+1]) {
294 cs -= fShellCrossSection.GetValueForComponent(Z, j, gammaEnergy);
295 }
296 if(cs <= 0.0 || j+1 == nn) { break; }
297 }
298 }
299 }
300
301 G4double bindingEnergy = (fParam[Z])[shellIdx*6 + 1];
302 //G4cout << "Z= " << Z << " shellIdx= " << shellIdx
303 // << " nShells= " << fNShells[Z]
304 // << " Ebind(keV)= " << bindingEnergy/keV
305 // << " Egamma(keV)= " << gammaEnergy/keV << G4endl;
306
307 const G4AtomicShell* shell = 0;
308
309 // no de-excitation from the last shell
310 if(fDeexcitationActive && shellIdx + 1 == nn) {
312 shell = fAtomDeexcitation->GetAtomicShell(Z, as);
313 }
314
315 // If binding energy of the selected shell is larger than photon energy
316 // do not generate secondaries
317 if(gammaEnergy < bindingEnergy) {
319 return;
320 }
321
322 // Primary outcoming electron
323 G4double eKineticEnergy = gammaEnergy - bindingEnergy;
324 G4double edep = bindingEnergy;
325
326 // Calculate direction of the photoelectron
327 G4ThreeVector electronDirection =
328 GetAngularDistribution()->SampleDirection(aDynamicGamma,
329 eKineticEnergy,
330 shellIdx,
331 couple->GetMaterial());
332
333 // The electron is created
334 G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
335 electronDirection,
336 eKineticEnergy);
337 fvect->push_back(electron);
338
339 // Sample deexcitation
340 if(shell) {
341 G4int index = couple->GetIndex();
342 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
343 size_t nbefore = fvect->size();
344
345 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
346 size_t nafter = fvect->size();
347 if(nafter > nbefore) {
348 for (size_t j=nbefore; j<nafter; ++j) {
349 edep -= ((*fvect)[j])->GetKineticEnergy();
350 }
351 }
352 }
353 }
354 // energy balance - excitation energy left
355 if(edep > 0.0) {
357 }
358}
G4AtomicShellEnumerator
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
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:131
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:508
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SetLimitNumberOfShells()

void G4LivermorePhotoElectricModel::SetLimitNumberOfShells ( G4int  n)
inline

Definition at line 110 of file G4LivermorePhotoElectricModel.hh.

111{
112 nShellLimit = n;
113}

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePhotoElectricModel::fParticleChange
protected

Definition at line 80 of file G4LivermorePhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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