Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4eBremsstrahlungRelModel.hh
이 파일의 문서화 페이지로 가기
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 // $Id: G4eBremsstrahlungRelModel.hh 108737 2018-03-02 13:49:56Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eBremsstrahlungRelModel
34 // extention of standard G4eBremsstrahlungModel
35 //
36 // Author: Andreas Schaelicke
37 //
38 // Creation date: 28.03.2008
39 //
40 // Modifications:
41 //
42 //
43 // Class Description:
44 //
45 // Implementation of energy loss for gamma emission by electrons and
46 // positrons including an improved version of the LPM effect
47 
48 // -------------------------------------------------------------------
49 //
50 
51 #ifndef G4eBremsstrahlungRelModel_h
52 #define G4eBremsstrahlungRelModel_h 1
53 
54 #include "G4VEmModel.hh"
55 #include "G4NistManager.hh"
56 #include "G4Exp.hh"
57 #include "G4Log.hh"
58 
60 class G4PhysicsVector;
61 
63 {
64 
65 public:
66 
68  const G4String& nam = "eBremLPM");
69 
70  ~G4eBremsstrahlungRelModel() = default;
71 
72  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
73 
74  virtual void InitialiseLocal(const G4ParticleDefinition*,
75  G4VEmModel* masterModel) override;
76 
78  const G4ParticleDefinition*,
80  G4double cutEnergy) override;
81 
83  G4double tkin,
85  G4double cutEnergy,
86  G4double maxEnergy = DBL_MAX) override;
87 
88  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89  const G4MaterialCutsCouple*,
90  const G4DynamicParticle*,
91  G4double cutEnergy,
92  G4double maxEnergy) override;
93 
94  virtual void SetupForMaterial(const G4ParticleDefinition*,
95  const G4Material*,G4double) override;
96 
97  virtual G4double MinPrimaryEnergy(const G4Material*,
98  const G4ParticleDefinition*,
99  G4double cut) override;
100 
101  inline void SetLPMconstant(G4double val);
102  inline G4double LPMconstant() const;
103 
104  inline void SetLowestKinEnergy(G4double);
105  inline G4double LowestKinEnergy() const;
106 
107 
108 protected:
109 
110  virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
111 
112  void SetCurrentElement(G4int);
113 
114 private:
115 
116  void InitialiseConstants();
117 
118  void CalcLPMFunctions(G4double gammaEnergy);
119 
120  G4double ComputeBremLoss(G4double cutEnergy);
121 
123 
125 
126  void SetParticle(const G4ParticleDefinition* p);
127 
128  inline G4double Phi1(G4double,G4double);
130  inline G4double Psi1(G4double,G4double);
132 
133  // hide assignment operator
136 
137 protected:
138 
143 
145 
146  // cash
152 
153  // scattering off electrons
156 
159 
161 
162 private:
163 
164  static const G4double xgi[8], wgi[8];
165  static const G4double Fel_light[5];
166  static const G4double Finel_light[5];
167 
168  // consts
175 
176  // cash
179 
180  // LPM effect
183 
184  // critical gamma energies
186 
187  // flags
189 
190 };
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
194 
196 {
197  // Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
198  return 20.863 - 2.*G4Log(1. + sqr(0.55846*gg) )
199  - 4.*( 1. - 0.6*G4Exp(-0.9*gg) - 0.4*G4Exp(-1.5*gg) );
200 }
201 
203 {
204  // Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
205  // return Phi1(gg,Z) -
206  return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
207 }
208 
210 {
211  // Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
212  return 28.340 - 2.*G4Log(1. + sqr(3.621*eps) )
213  - 4.*( 1. - 0.7*G4Exp(-8*eps) - 0.3*G4Exp(-29.*eps) );
214 }
215 
217 {
218  // Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
219  return 2./(3.*(1. + 40.*eps +400.*eps*eps) );
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
223 
224 inline
226 {
227  fLPMconstant = val;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
232 inline
234 {
235  return fLPMconstant;
236 }
237 
239 {
240  lowestKinEnergy = val;
241 }
242 
243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244 
246 {
247  return lowestKinEnergy;
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
251 
252 
253 #endif
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
const G4ParticleDefinition * particle
G4double ComputeXSectionPerAtom(G4double cutEnergy)
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
G4double Psi1(G4double, G4double)
const char * p
Definition: xmltok.h:285
G4double Phi1(G4double, G4double)
static const G4double Finel_light[5]
G4double ComputeBremLoss(G4double cutEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
G4double Psi1M2(G4double, G4double)
~G4eBremsstrahlungRelModel()=default
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4eBremsstrahlungRelModel & operator=(const G4eBremsstrahlungRelModel &right)=delete
void CalcLPMFunctions(G4double gammaEnergy)
int G4int
Definition: G4Types.hh:78
void SetParticle(const G4ParticleDefinition *p)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static const G4double Fel_light[5]
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
T sqr(const T &x)
Definition: templates.hh:145
G4ParticleChangeForLoss * fParticleChange
G4double Phi1M2(G4double, G4double)
static const G4double eps
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut) override
#define DBL_MAX
Definition: templates.hh:83
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override