Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eBremsstrahlungRelModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eBremsstrahlungRelModel
33// extention of standard G4eBremsstrahlungModel
34//
35// Author: Andreas Schaelicke
36//
37// Creation date: 28.03.2008
38//
39// Modifications:
40//
41// 15.07.18 introduced data structures to store LPM functions and element depen
42// dent data for faster run-time computation (see more in .cc M.Novak)
43//
44// Class Description:
45//
46// Implementation of energy loss for gamma emission by electrons and
47// positrons including an improved version of the LPM effect
48
49// -------------------------------------------------------------------
50//
51
52#ifndef G4eBremsstrahlungRelModel_h
53#define G4eBremsstrahlungRelModel_h 1
54
55#include "G4VEmModel.hh"
56
58
60
61public:
62
63 explicit G4eBremsstrahlungRelModel(const G4ParticleDefinition* p=nullptr,
64 const G4String& nam="eBremLPM");
65
67
68 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
69
71 G4VEmModel* masterModel) override;
72
75 G4double ekin,
76 G4double cutEnergy) override;
77
79 G4double ekin,
80 G4double zet,
82 G4double cutEnergy,
83 G4double maxEnergy = DBL_MAX) override;
84
85 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
87 const G4DynamicParticle*,
88 G4double cutEnergy,
89 G4double maxEnergy) override;
90
92 const G4Material*, G4double) override;
93
96 G4double cutEnergy) override;
97
98protected:
99
100 virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
101
102 void SetParticle(const G4ParticleDefinition* p);
103
104private:
105
106 G4double ComputeBremLoss(G4double cutEnergy);
107
108 G4double ComputeXSectionPerAtom(G4double cutEnergy);
109
110 G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
111
112 // init special data per element i.e. per Z
113 void InitialiseElementData();
114
115 // methods for initialisation and run-time evaluation of LPM functions:
116 void InitLPMFunctions();
117
118 void ComputeLPMfunctions(G4double& funcXiS,
119 G4double& funcGS,
120 G4double& funcPhiS,
121 const G4double egamma);
122
123 void GetLPMFunctions(G4double& lpmGs,
124 G4double& lpmPhis,
125 const G4double ss);
126
127 void ComputeLPMGsPhis(G4double& funcGS,
128 G4double& funcPhiS,
129 const G4double varShat);
130
131 // for evaluating screening related functions
132 void ComputeScreeningFunctions(G4double& phi1,
133 G4double& phi1m2,
134 G4double& psi1,
135 G4double& psi1m2,
136 const G4double gam,
137 const G4double eps);
138 // hide assignment operator and cctr
140 (const G4eBremsstrahlungRelModel& right) = delete;
142
143private:
144
145 G4bool fIsUseCompleteScreening = false;
146 G4bool fIsInitializer = false;
147 G4bool fUseLPM = true;
148
149protected:
150
154 //
159 // cash
166 // scattering off electrons
169
170private:
171
172 // LPM related members
173 G4double fLPMEnergyThreshold;
174 G4double fLPMEnergy;
175
176protected:
177
178 static const G4double gBremFactor;
180
181private:
182
183 static const G4int gMaxZet;
184 //
185 static const G4double gLPMconstant;
186 //
187 static const G4double gXGL[8];
188 static const G4double gWGL[8];
189 static const G4double gFelLowZet[8];
190 static const G4double gFinelLowZet[8];
191 //
192 struct ElementData {
193 /** @brief \f$ \ln(Z) \f$ */
194 G4double fLogZ;
195 /** @brief \f$ \ln(Z)/3 + f_c \f$ */
196 G4double fFz;
197 /** @brief \f$ ((Fel-fc)+Finel*invZ)\f$ */
198 G4double fZFactor1;
199 /** @brief \f$ (Fel-fc)\f$ */
200 G4double fZFactor11;
201 /** @brief \f$ (1.0+invZ)/12 \f$ */
202 G4double fZFactor2;
203 // LPM variables
204 G4double fVarS1;
205 G4double fILVarS1;
206 G4double fILVarS1Cond;
207 // constant factors to the screening function evaluations
208 G4double fGammaFactor;
209 G4double fEpsilonFactor;
210 };
211 //
212 struct LPMFuncs {
213 LPMFuncs() : fIsInitialized(false), fISDelta(100.), fSLimit(2.) {}
214 G4bool fIsInitialized;
215 G4double fISDelta;
216 G4double fSLimit;
217 std::vector<G4double> fLPMFuncG;
218 std::vector<G4double> fLPMFuncPhi;
219 };
220 //
221 static LPMFuncs gLPMFuncs;
222 static std::vector<ElementData*> gElementData;
223
224};
225
226#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetParticle(const G4ParticleDefinition *p)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
const G4ParticleDefinition * fPrimaryParticle
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremLPM")
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4ParticleChangeForLoss * fParticleChange
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
#define DBL_MAX
Definition templates.hh:62