Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
64public:
65
67 const G4String& nam = "BetheHeitlerLPM");
68
70
71 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
72
75 G4double kinEnergy,
76 G4double Z,
77 G4double A=0.,
78 G4double cut=0.,
79 G4double emax=DBL_MAX);
80
81 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
83 const G4DynamicParticle*,
84 G4double tmin,
85 G4double maxEnergy);
86
87 virtual void SetupForMaterial(const G4ParticleDefinition*,
88 const G4Material*,G4double);
89
90 // * fast inline functions *
91 inline void SetCurrentElement(G4double /*Z*/);
92
93 // set / get methods
94 inline void SetLPMconstant(G4double val);
95 inline G4double LPMconstant() const;
96
97 inline void SetLPMflag(G4bool);
98 inline G4bool LPMflag() const;
99
100protected:
101 // screening functions
102 inline G4double Phi1(G4double delta) const;
103 inline G4double Phi2(G4double delta) const;
104 inline G4double DeltaMax() const;
105 inline G4double DeltaMin(G4double) const;
106 // lpm functions
107 void CalcLPMFunctions(G4double k, G4double eplus);
108
109 // obsolete
110 G4double ScreenFunction1(G4double ScreenVariable);
111 G4double ScreenFunction2(G4double ScreenVariable);
112
114
115 G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
116 G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
117
118 // hide assignment operator
121
123
128
131
132 // cash
136
137 // LPM effect
140
141 // consts
143
144 static const G4double xgi[8], wgi[8];
145 static const G4double Fel_light[5];
146 static const G4double Finel_light[5];
147 static const G4double facFel;
148 static const G4double facFinel;
149
150 static const G4double preS1, logTwo;
151
152};
153
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
157inline
159{
160 fLPMconstant = val;
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
165inline
167{
168 return fLPMconstant;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
173inline
175{
176 fLPMflag = val;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181inline
183{
184 return fLPMflag;
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188
190{
191 if(Z != currentZ) {
192 currentZ = Z;
193
194 G4int iz = G4int(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 Fel = facFel - lnZ/3. ;
205 Finel = facFinel - 2.*lnZ/3. ;
206 }
208 }
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
212
214{
215 G4double screenVal;
216
217 if (delta > 1.)
218 screenVal = 21.12 - 4.184*std::log(delta+0.952);
219 else
220 screenVal = 20.868 - delta*(3.242 - 0.625*delta);
221
222 return screenVal;
223}
224
225//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
226
228{
229 G4double screenVal;
230
231 if (delta > 1.)
232 screenVal = 21.12 - 4.184*std::log(delta+0.952);
233 else
234 screenVal = 20.209 - delta*(1.930 + 0.086*delta);
235
236 return screenVal;
237}
238
240
241// compute the value of the screening function 3*PHI1 - PHI2
242
243{
244 G4double screenVal;
245
246 if (ScreenVariable > 1.)
247 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
248 else
249 screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
250
251 return screenVal;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255
257
258// compute the value of the screening function 1.5*PHI1 + 0.5*PHI2
259
260{
261 G4double screenVal;
262
263 if (ScreenVariable > 1.)
264 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
265 else
266 screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
267
268 return screenVal;
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272
273
275{
276 // k > 50 MeV
277 G4double FZ = 8.*(lnZ/3. + fCoulomb);
278 return std::exp( (42.24-FZ)/8.368 ) + 0.952;
279}
280
282{
283 return 4.*136./z13*(CLHEP::electron_mass_c2/k);
284}
285
286
287
288#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4double GetZ13(G4double Z)
G4double GetLOGZ(G4int Z)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleDefinition * thePositron
static const G4double Finel_light[5]
static const G4double wgi[8]
G4double ScreenFunction2(G4double ScreenVariable)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
G4ParticleDefinition * theElectron
static const G4double Fel_light[5]
static const G4double xgi[8]
G4PairProductionRelModel & operator=(const G4PairProductionRelModel &right)
G4PairProductionRelModel(const G4PairProductionRelModel &)
void CalcLPMFunctions(G4double k, G4double eplus)
G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)
G4double ScreenFunction1(G4double ScreenVariable)
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double Phi2(G4double delta) const
G4ParticleChangeForGamma * fParticleChange
G4double DeltaMin(G4double) const
G4double Phi1(G4double delta) const
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:391
#define DBL_MAX
Definition: templates.hh:83