Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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$
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 
58 class G4PhysicsVector;
59 
61 {
62 
63 public:
64 
66  const G4String& nam = "eBremLPM");
67 
69 
70  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
71 
73  const G4ParticleDefinition*,
74  G4double kineticEnergy,
75  G4double cutEnergy);
76 
78  G4double tkin,
80  G4double cutEnergy,
81  G4double maxEnergy = DBL_MAX);
82 
83  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
84  const G4MaterialCutsCouple*,
85  const G4DynamicParticle*,
86  G4double cutEnergy,
87  G4double maxEnergy);
88 
89  virtual void SetupForMaterial(const G4ParticleDefinition*,
90  const G4Material*,G4double);
91 
92  inline void SetLPMconstant(G4double val);
93  inline G4double LPMconstant() const;
94 
95 protected:
96 
97  virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
98 
99  // * fast inline functions *
100  inline void SetCurrentElement(const G4double);
101 
102 private:
103 
104  void InitialiseConstants();
105 
106  void CalcLPMFunctions(G4double gammaEnergy);
107 
108  G4double ComputeBremLoss(G4double cutEnergy);
109 
110  G4double ComputeXSectionPerAtom(G4double cutEnergy);
111 
112  G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
113 
114  void SetParticle(const G4ParticleDefinition* p);
115 
116  inline G4double Phi1(G4double,G4double);
117  inline G4double Phi1M2(G4double,G4double);
118  inline G4double Psi1(G4double,G4double);
119  inline G4double Psi1M2(G4double,G4double);
120 
121  // hide assignment operator
124 
125 protected:
126 
131 
133 
134  // cash
139  //G4double z13, z23, lnZ;
140  //G4double Fel, Finel, fCoulomb, fMax;
143 
145 
146 private:
147 
148  static const G4double xgi[8], wgi[8];
149  static const G4double Fel_light[5];
150  static const G4double Finel_light[5];
151 
152  // consts
153  G4double lowKinEnergy;
154  G4double fMigdalConstant;
155  G4double fLPMconstant;
156  G4double energyThresholdLPM;
157  G4double facFel, facFinel;
158  G4double preS1,logTwo;
159 
160  // cash
161  //G4double particleMass;
162  //G4double kinEnergy;
163  //G4double totalEnergy;
164  //G4double currentZ;
165  G4double z13, z23, lnZ;
166  G4double Fel, Finel, fCoulomb, fMax;
167  //G4double densityFactor;
168  //G4double densityCorr;
169 
170  // LPM effect
171  G4double lpmEnergy;
172  G4PhysicsVector *fXiLPM, *fPhiLPM, *fGLPM;
173  G4double xiLPM, phiLPM, gLPM;
174 
175  // critical gamma energies
176  G4double klpm, kp;
177 
178  // flags
179  G4bool use_completescreening;
180  G4bool isInitialised;
181 };
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
186 {
187  if(Z != currentZ) {
188  currentZ = Z;
189 
190  G4int iz = G4int(Z);
191  z13 = nist->GetZ13(iz);
192  z23 = z13*z13;
193  lnZ = nist->GetLOGZ(iz);
194 
195  if (iz <= 4) {
196  Fel = Fel_light[iz];
197  Finel = Finel_light[iz] ;
198  }
199  else {
200  Fel = facFel - lnZ/3. ;
201  Finel = facFinel - 2.*lnZ/3. ;
202  }
203 
204  fCoulomb = GetCurrentElement()->GetfCoulomb();
205  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
206  }
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210 
211 
212 inline G4double G4eBremsstrahlungRelModel::Phi1(G4double gg, G4double)
213 {
214  // Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
215  return 20.863 - 2.*std::log(1. + sqr(0.55846*gg) )
216  - 4.*( 1. - 0.6*std::exp(-0.9*gg) - 0.4*std::exp(-1.5*gg) );
217 }
218 
219 inline G4double G4eBremsstrahlungRelModel::Phi1M2(G4double gg, G4double)
220 {
221  // Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
222  // return Phi1(gg,Z) -
223  return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
224 }
225 
226 inline G4double G4eBremsstrahlungRelModel::Psi1(G4double eps, G4double)
227 {
228  // Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
229  return 28.340 - 2.*std::log(1. + sqr(3.621*eps) )
230  - 4.*( 1. - 0.7*std::exp(-8*eps) - 0.3*std::exp(-29.*eps) );
231 }
232 
233 inline G4double G4eBremsstrahlungRelModel::Psi1M2(G4double eps, G4double)
234 {
235  // Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
236  return 2./(3.*(1. + 40.*eps +400.*eps*eps) );
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
240 
241 inline
243 {
244  fLPMconstant = val;
245 }
246 
247 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
248 
249 inline
251 {
252  return fLPMconstant;
253 }
254 
255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256 
257 
258 #endif