Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PairProductionRelModel.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: G4PairProductionRelModel.hh 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4PairProductionRelModel
34 //
35 // Author: Andreas Schaelicke
36 //
37 // Creation date: 02.04.2009
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Implementation of gamma convertion to e+e- in the field of a nucleus
44 // relativistic approximation
45 //
46 
47 // -------------------------------------------------------------------
48 //
49 
50 #ifndef G4PairProductionRelModel_h
51 #define G4PairProductionRelModel_h 1
52 
54 
55 #include "G4VEmModel.hh"
56 #include "G4PhysicsTable.hh"
57 #include "G4NistManager.hh"
58 
60 
62 {
63 
64 public:
65 
66  explicit G4PairProductionRelModel(const G4ParticleDefinition* p = 0,
67  const G4String& nam = "BetheHeitlerLPM");
68 
69  virtual ~G4PairProductionRelModel();
70 
71  virtual void Initialise(const G4ParticleDefinition*,
72  const G4DataVector&) override;
73 
74  virtual void InitialiseLocal(const G4ParticleDefinition*,
75  G4VEmModel* masterModel) override;
76 
78  const G4ParticleDefinition*,
79  G4double kinEnergy,
80  G4double Z,
81  G4double A=0.,
82  G4double cut=0.,
83  G4double emax=DBL_MAX) override;
84 
85  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
86  const G4MaterialCutsCouple*,
87  const G4DynamicParticle*,
88  G4double tmin,
89  G4double maxEnergy) override;
90 
91  virtual void SetupForMaterial(const G4ParticleDefinition*,
92  const G4Material*,G4double) override;
93 
94  // * fast inline functions *
95  inline void SetCurrentElement(G4double Z);
96 
97  // set / get methods
98  inline void SetLPMconstant(G4double val);
99  inline G4double LPMconstant() const;
100 
101  inline void SetLPMflag(G4bool);
102  inline G4bool LPMflag() const;
103 
104 protected:
105 
106  // screening functions
107  inline G4double Phi1(G4double delta) const;
108  inline G4double Phi2(G4double delta) const;
109  inline G4double ScreenFunction1(G4double ScreenVariable);
110  inline G4double ScreenFunction2(G4double ScreenVariable);
111  inline G4double DeltaMax() const;
112  inline G4double DeltaMin(G4double) const;
113 
114  // lpm functions
116 
118 
119  G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
120  G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
121 
122  // hide assignment operator
123  G4PairProductionRelModel & operator=
124  (const G4PairProductionRelModel &right) = delete;
126 
128 
133 
136 
137  // cash
141 
142  // LPM effect
145 
146  // consts
148 
149  static const G4double xgi[8], wgi[8];
150  static const G4double Fel_light[5];
151  static const G4double Finel_light[5];
152  static const G4double facFel;
153  static const G4double facFinel;
154 
156 
157 };
158 
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161 
162 inline
164 {
165  fLPMconstant = val;
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
170 inline
172 {
173  return fLPMconstant;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177 
178 inline
180 {
181  fLPMflag = val;
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185 
186 inline
188 {
189  return fLPMflag;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
195 {
196  if(Z != currentZ) {
197  currentZ = Z;
198 
199  G4int iz = G4lrint(Z);
200  z13 = nist->GetZ13(iz);
201  z23 = z13*z13;
202  lnZ = nist->GetLOGZ(iz);
203 
204  if (iz <= 4) {
205  Fel = Fel_light[iz];
206  Finel = Finel_light[iz] ;
207  }
208  else {
209  Fel = facFel - lnZ/3. ;
210  Finel = facFinel - 2.*lnZ/3. ;
211  }
213  }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217 
219 {
220  G4double screenVal;
221 
222  if (delta > 1.)
223  screenVal = 21.12 - 4.184*std::log(delta+0.952);
224  else
225  screenVal = 20.868 - delta*(3.242 - 0.625*delta);
226 
227  return screenVal;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
233 {
234  G4double screenVal;
235 
236  if (delta > 1.)
237  screenVal = 21.12 - 4.184*std::log(delta+0.952);
238  else
239  screenVal = 20.209 - delta*(1.930 + 0.086*delta);
240 
241  return screenVal;
242 }
243 
245 
246 // compute the value of the screening function 3*PHI1 - PHI2
247 
248 {
249  G4double screenVal;
250 
251  if (ScreenVariable > 1.)
252  screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
253  else
254  screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
255 
256  return screenVal;
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260 
262 
263 // compute the value of the screening function 1.5*PHI1 + 0.5*PHI2
264 
265 {
266  G4double screenVal;
267 
268  if (ScreenVariable > 1.)
269  screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
270  else
271  screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
272 
273  return screenVal;
274 }
275 
276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
277 
278 
280 {
281  // k > 50 MeV
282  G4double FZ = 8.*(lnZ/3. + fCoulomb);
283  return std::exp( (42.24-FZ)/8.368 ) + 0.952;
284 }
285 
287 {
288  return 4.*136./z13*(CLHEP::electron_mass_c2/k);
289 }
290 
291 
292 
293 #endif
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4double Phi1(G4double delta) const
G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
const char * p
Definition: xmltok.h:285
G4double Phi2(G4double delta) const
G4PairProductionRelModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitlerLPM")
G4double GetfCoulomb() const
Definition: G4Element.hh:191
static const G4double Finel_light[5]
G4ParticleDefinition * thePositron
int G4int
Definition: G4Types.hh:78
static constexpr double electron_mass_c2
static const G4double xgi[8]
G4double ScreenFunction1(G4double ScreenVariable)
G4double GetZ13(G4double Z) const
double A(double temperature)
static const G4double wgi[8]
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
Definition: G4SIunits.hh:199
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double DeltaMin(G4double) const
void CalcLPMFunctions(G4double k, G4double eplus)
G4double GetLOGZ(G4int Z) const
G4double ScreenFunction2(G4double ScreenVariable)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
static const G4double emax
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4ParticleChangeForGamma * fParticleChange
int G4lrint(double ad)
Definition: templates.hh:163
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const G4double Fel_light[5]
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466