Geant4 10.7.0
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
64 const G4String& nam="eBremLPM");
65
67
68 virtual void Initialise(const G4ParticleDefinition*,
69 const G4DataVector&) override;
70
71 virtual void InitialiseLocal(const G4ParticleDefinition*,
72 G4VEmModel* masterModel) override;
73
76 G4double ekin,
77 G4double cutEnergy) override;
78
80 G4double ekin,
81 G4double zet,
83 G4double cutEnergy,
84 G4double maxEnergy = DBL_MAX) override;
85
86 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
88 const G4DynamicParticle*,
89 G4double cutEnergy,
90 G4double maxEnergy) override;
91
92 virtual void SetupForMaterial(const G4ParticleDefinition*,
93 const G4Material*,
94 G4double) override;
95
98 G4double cutEnergy) override;
99
100protected:
101
102 virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
103
104 void SetParticle(const G4ParticleDefinition* p);
105
106private:
107
108 G4double ComputeBremLoss(G4double cutEnergy);
109
110 G4double ComputeXSectionPerAtom(G4double cutEnergy);
111
112 G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
113
114 // init special data per element i.e. per Z
115 void InitialiseElementData();
116
117 // methods for initialisation and run-time evaluation of LPM functions:
118 void InitLPMFunctions();
119
120 void ComputeLPMfunctions(G4double& funcXiS,
121 G4double& funcGS,
122 G4double& funcPhiS,
123 const G4double egamma);
124
125 void GetLPMFunctions(G4double& lpmGs,
126 G4double& lpmPhis,
127 const G4double s);
128
129 void ComputeLPMGsPhis(G4double& funcGS,
130 G4double& funcPhiS,
131 const G4double varShat);
132 //
133 // for evaluating screening related functions
134 void ComputeScreeningFunctions(G4double& phi1,
135 G4double& phi1m2,
136 G4double& psi1,
137 G4double& psi1m2,
138 const G4double gam,
139 const G4double eps);
140 // hide assignment operator and cctr
141 G4eBremsstrahlungRelModel& operator=(const G4eBremsstrahlungRelModel& right) = delete;
143
144protected:
145
149 //
151 // cash
158 // scattering off electrons
161 //
162 static const G4double gBremFactor;
164 //
168
169private:
170 static const G4int gMaxZet;
171 //
172 static const G4double gLPMconstant;
173 //
174 static const G4double gXGL[8];
175 static const G4double gWGL[8];
176 static const G4double gFelLowZet[8];
177 static const G4double gFinelLowZet[8];
178 //
179 struct ElementData {
180 /** @brief \f$ \ln(Z) \f$ */
181 G4double fLogZ;
182 /** @brief \f$ \ln(Z)/3 + f_c \f$ */
183 G4double fFz;
184 /** @brief \f$ ((Fel-fc)+Finel*invZ)\f$ */
185 G4double fZFactor1;
186 /** @brief \f$ (Fel-fc)\f$ */
187 G4double fZFactor11;
188 /** @brief \f$ (1.0+invZ)/12 \f$ */
189 G4double fZFactor2;
190 // LPM variables
191 G4double fVarS1;
192 G4double fILVarS1;
193 G4double fILVarS1Cond;
194 // constant factors to the screening function evaluations
195 G4double fGammaFactor;
196 G4double fEpsilonFactor;
197 };
198 //
199 struct LPMFuncs {
200 LPMFuncs() : fIsInitialized(false), fISDelta(100.), fSLimit(2.) {}
201 G4bool fIsInitialized;
202 G4double fISDelta;
203 G4double fSLimit;
204 std::vector<G4double> fLPMFuncG;
205 std::vector<G4double> fLPMFuncPhi;
206 };
207 //
208 static LPMFuncs gLPMFuncs;
209 static std::vector<ElementData*> gElementData;
210 //
211 G4bool fIsUseCompleteScreening;
212 // LPM related members
213 G4double fLPMEnergyThreshold;
214 G4double fLPMEnergy;
215
216};
217
218#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetParticle(const G4ParticleDefinition *p)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
const G4ParticleDefinition * fPrimaryParticle
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
G4ParticleDefinition * fGammaParticle
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4ParticleChangeForLoss * fParticleChange
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
#define DBL_MAX
Definition: templates.hh:62