Geant4  10.02.p02
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 // $Id: G4eBremsstrahlungRelModel.hh 93567 2015-10-26 14:51:41Z 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 
71 
72  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
73 
74  virtual void InitialiseLocal(const G4ParticleDefinition*,
75  G4VEmModel* masterModel);
76 
78  const G4ParticleDefinition*,
79  G4double kineticEnergy,
80  G4double cutEnergy);
81 
83  G4double tkin,
84  G4double Z, G4double,
85  G4double cutEnergy,
86  G4double maxEnergy = DBL_MAX);
87 
88  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89  const G4MaterialCutsCouple*,
90  const G4DynamicParticle*,
91  G4double cutEnergy,
92  G4double maxEnergy);
93 
94  virtual void SetupForMaterial(const G4ParticleDefinition*,
95  const G4Material*,G4double);
96 
97  virtual G4double MinPrimaryEnergy(const G4Material*,
98  const G4ParticleDefinition*,
99  G4double cut);
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  // * fast inline functions *
113  inline void SetCurrentElement(const G4double);
114 
115 private:
116 
117  void InitialiseConstants();
118 
119  void CalcLPMFunctions(G4double gammaEnergy);
120 
121  G4double ComputeBremLoss(G4double cutEnergy);
122 
124 
126 
127  void SetParticle(const G4ParticleDefinition* p);
128 
129  inline G4double Phi1(G4double,G4double);
131  inline G4double Psi1(G4double,G4double);
133 
134  // hide assignment operator
137 
138 protected:
139 
144 
146 
147  // cash
154 
156 
157 private:
158 
159  static const G4double xgi[8], wgi[8];
160  static const G4double Fel_light[5];
161  static const G4double Finel_light[5];
162 
163  // consts
170 
171  // cash
174 
175  // LPM effect
179 
180  // critical gamma energies
182 
183  // flags
185 };
186 
187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188 
190 {
191  if(Z != currentZ) {
192  currentZ = Z;
193 
194  G4int iz = G4lrint(Z);
195  z13 = nist->GetZ13(iz);
196  z23 = z13*z13;
197  lnZ = nist->GetLOGZ(iz);
198 
199  if (iz <= 4) {
200  Fel = Fel_light[iz];
201  Finel = Finel_light[iz] ;
202  }
203  else {
204  G4double lnzt = lnZ/3.;
205  Fel = facFel - lnzt;
206  Finel = facFinel - 2*lnzt;
207  }
208 
210  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
211  }
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215 
216 
218 {
219  // Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
220  return 20.863 - 2.*G4Log(1. + sqr(0.55846*gg) )
221  - 4.*( 1. - 0.6*G4Exp(-0.9*gg) - 0.4*G4Exp(-1.5*gg) );
222 }
223 
225 {
226  // Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
227  // return Phi1(gg,Z) -
228  return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
229 }
230 
232 {
233  // Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
234  return 28.340 - 2.*G4Log(1. + sqr(3.621*eps) )
235  - 4.*( 1. - 0.7*G4Exp(-8*eps) - 0.3*G4Exp(-29.*eps) );
236 }
237 
239 {
240  // Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
241  return 2./(3.*(1. + 40.*eps +400.*eps*eps) );
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245 
246 inline
248 {
249  fLPMconstant = val;
250 }
251 
252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253 
254 inline
256 {
257  return fLPMconstant;
258 }
259 
261 {
262  lowestKinEnergy = val;
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266 
268 {
269  return lowestKinEnergy;
270 }
271 
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
273 
274 
275 #endif
G4double Psi1(G4double, G4double)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double GetZ13(G4double Z)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
static const G4double eps
G4eBremsstrahlungRelModel & operator=(const G4eBremsstrahlungRelModel &right)
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double Psi1M2(G4double, G4double)
int G4int
Definition: G4Types.hh:78
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
static const G4double Finel_light[5]
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
G4double Phi1(G4double, G4double)
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut)
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double Fel_light[5]
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
G4double Phi1M2(G4double, G4double)
G4double ComputeXSectionPerAtom(G4double cutEnergy)
G4double GetLOGZ(G4int Z)
G4ParticleChangeForLoss * fParticleChange
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetParticle(const G4ParticleDefinition *p)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4double ComputeBremLoss(G4double cutEnergy)
void CalcLPMFunctions(G4double gammaEnergy)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:467