source: trunk/source/processes/electromagnetic/standard/include/G4eBremsstrahlungRelModel.hh @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 7.3 KB
Line 
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,v 1.14 2010/10/26 10:35:22 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4eBremsstrahlungRelModel
35//                extention of standard G4eBremsstrahlungModel
36//
37// Author:        Andreas Schaelicke
38//
39// Creation date: 28.03.2008
40//
41// Modifications:
42//
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#include "G4NistManager.hh"
57
58class G4ParticleChangeForLoss;
59class G4PhysicsVector;
60
61class G4eBremsstrahlungRelModel : public G4VEmModel
62{
63
64public:
65
66  G4eBremsstrahlungRelModel(const G4ParticleDefinition* p = 0, 
67                            const G4String& nam = "eBremLPM");
68
69  virtual ~G4eBremsstrahlungRelModel();
70
71  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
72
73  virtual G4double MinEnergyCut(const G4ParticleDefinition*, 
74                                const G4MaterialCutsCouple*);
75
76  virtual G4double ComputeDEDXPerVolume(const G4Material*,
77                                        const G4ParticleDefinition*,
78                                        G4double kineticEnergy,
79                                        G4double cutEnergy);
80                                       
81  virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
82                                              G4double tkin, 
83                                              G4double Z,   G4double,
84                                              G4double cutEnergy,
85                                              G4double maxEnergy = DBL_MAX);
86 
87  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
88                                 const G4MaterialCutsCouple*,
89                                 const G4DynamicParticle*,
90                                 G4double cutEnergy,
91                                 G4double maxEnergy);
92
93  virtual void SetupForMaterial(const G4ParticleDefinition*,
94                                const G4Material*,G4double);
95
96  inline void SetLPMconstant(G4double val);
97  inline G4double LPMconstant() const;
98
99private:
100
101  void InitialiseConstants();
102
103  void CalcLPMFunctions(G4double gammaEnergy);
104
105  G4double ComputeBremLoss(G4double cutEnergy);
106
107  G4double ComputeXSectionPerAtom(G4double cutEnergy);
108
109  G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
110
111  G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
112
113  void SetParticle(const G4ParticleDefinition* p);
114
115  // * fast inline functions *
116  inline void SetCurrentElement(const G4double);
117
118  inline G4double Phi1(G4double,G4double);
119  inline G4double Phi1M2(G4double,G4double);
120  inline G4double Psi1(G4double,G4double);
121  inline G4double Psi1M2(G4double,G4double);
122
123  // hide assignment operator
124  G4eBremsstrahlungRelModel & operator=(const  G4eBremsstrahlungRelModel &right);
125  G4eBremsstrahlungRelModel(const  G4eBremsstrahlungRelModel&);
126
127protected:
128
129  G4NistManager*              nist;
130  const G4ParticleDefinition* particle;
131  G4ParticleDefinition*       theGamma;
132  G4ParticleChangeForLoss*    fParticleChange;
133
134  static const G4double xgi[8], wgi[8];
135  static const G4double Fel_light[5];
136  static const G4double Finel_light[5];
137
138  G4double minThreshold;
139
140  // cash
141  G4double particleMass;
142  G4double kinEnergy;
143  G4double totalEnergy;
144  G4double currentZ;
145  G4double z13, z23, lnZ;
146  G4double Fel, Finel, fCoulomb, fMax; 
147  G4double densityFactor;
148  G4double densityCorr;
149
150  // LPM effect
151  G4double lpmEnergy;
152  G4PhysicsVector  *fXiLPM, *fPhiLPM, *fGLPM;
153  G4double xiLPM, phiLPM, gLPM;
154
155  // critical gamma energies
156  G4double klpm, kp;
157
158  G4bool   isElectron;
159
160private:
161
162  // consts
163  G4double lowKinEnergy;
164  G4double fMigdalConstant;
165  G4double fLPMconstant;
166  G4double bremFactor;
167  G4double energyThresholdLPM;
168  G4double facFel, facFinel;
169  G4double preS1,logTwo;
170
171  G4bool   use_completescreening;
172
173  G4bool   isInitialised;
174};
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
178inline void G4eBremsstrahlungRelModel::SetCurrentElement(const G4double Z)
179{
180  if(Z != currentZ) {
181    currentZ = Z;
182
183    G4int iz = G4int(Z);
184    z13 = nist->GetZ13(iz);
185    z23 = z13*z13;
186    lnZ = nist->GetLOGZ(iz);
187
188    if (iz <= 4) {
189      Fel = Fel_light[iz]; 
190      Finel = Finel_light[iz] ; 
191    }
192    else {
193      Fel = facFel - lnZ/3. ;
194      Finel = facFinel - 2.*lnZ/3. ;
195    }
196
197    fCoulomb = GetCurrentElement()->GetfCoulomb();
198    fMax = Fel-fCoulomb + Finel/currentZ  +  (1.+1./currentZ)/12.;
199  }
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203
204
205inline G4double G4eBremsstrahlungRelModel::Phi1(G4double gg, G4double)
206{
207  //       Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
208  return 20.863 - 2.*std::log(1. + sqr(0.55846*gg) )
209    - 4.*( 1. - 0.6*std::exp(-0.9*gg) - 0.4*std::exp(-1.5*gg) );
210}
211
212inline G4double G4eBremsstrahlungRelModel::Phi1M2(G4double gg, G4double)
213{
214  //       Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
215  // return Phi1(gg,Z) -
216  return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
217}
218
219inline G4double G4eBremsstrahlungRelModel::Psi1(G4double eps, G4double)
220{
221  //       Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
222  return 28.340 - 2.*std::log(1. + sqr(3.621*eps) )
223    - 4.*( 1. - 0.7*std::exp(-8*eps) - 0.3*std::exp(-29.*eps) );
224}
225
226inline G4double G4eBremsstrahlungRelModel::Psi1M2(G4double eps, G4double)
227{
228  //       Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
229  return  2./(3.*(1. + 40.*eps +400.*eps*eps) );
230}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
233
234inline 
235void G4eBremsstrahlungRelModel::SetLPMconstant(G4double val) 
236{
237  fLPMconstant = val;
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241
242inline 
243G4double G4eBremsstrahlungRelModel::LPMconstant() const 
244{
245  return fLPMconstant;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
250
251#endif
Note: See TracBrowser for help on using the repository browser.