Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PairProductionRelModel.hh
이 파일의 문서화 페이지로 가기
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 110527 2018-05-29 06:09:58Z 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 // 28-05-18 New version with improved screening function approximation, improved
41 // LPM function approximation, efficiency, documentation and cleanup.
42 // Corrected call to selecting target atom in the final state sampling.
43 // (M. Novak)
44 //
45 // Class Description:
46 //
47 // Implementation of gamma convertion to e+e- in the field of a nucleus
48 // relativistic approximation
49 //
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4PairProductionRelModel_h
55 #define G4PairProductionRelModel_h 1
56 
58 
59 #include "G4VEmModel.hh"
60 #include "G4Log.hh"
61 #include "G4Exp.hh"
62 #include "G4Pow.hh"
63 
64 #include <vector>
65 
67 
69 {
70 
71 public:
72 
73  explicit G4PairProductionRelModel(const G4ParticleDefinition* p = nullptr,
74  const G4String& nam = "BetheHeitlerLPM");
75 
76  virtual ~G4PairProductionRelModel();
77 
78  virtual void Initialise(const G4ParticleDefinition*,
79  const G4DataVector&) override;
80 
81  virtual void InitialiseLocal(const G4ParticleDefinition*,
82  G4VEmModel* masterModel) override;
83 
85  const G4ParticleDefinition*,
86  G4double kinEnergy,
87  G4double Z,
88  G4double A=0.,
89  G4double cut=0.,
90  G4double emax=DBL_MAX) override;
91 
92  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
93  const G4MaterialCutsCouple*,
94  const G4DynamicParticle*,
95  G4double tmin,
96  G4double maxEnergy) override;
97 
98  virtual void SetupForMaterial(const G4ParticleDefinition*,
99  const G4Material*,G4double) override;
100 
101  inline void SetLPMflag(G4bool val) { fIsUseLPMCorrection = val; }
102  inline G4bool LPMflag() const { return fIsUseLPMCorrection; }
103 
104 protected:
105 
106  // for evaluating screening related functions
107  inline void ComputePhi12(const G4double delta, G4double &phi1, G4double &phi2);
108  inline G4double ScreenFunction1(const G4double delta);
109  inline G4double ScreenFunction2(const G4double delta);
110  inline void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2);
111  // helper methods for cross-section computation under different approximations
113  G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy,
114  G4double Z);
115  G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy,
116  G4double Z);
117 
118 
119 private:
120 
121  // for creating some data structure per Z with often used comp. intensive data
122  void InitialiseElementData();
123  struct ElementData {
132  };
133  // for precomputing comp. intensive parts of LPM suppression functions and
134  // using them at run-time
135  void InitLPMFunctions();
136  void ComputeLPMGsPhis(G4double &funcGS, G4double &funcPhiS,
137  const G4double varShat);
138  void GetLPMFunctions(G4double &lpmGs, G4double &lpmPhis, const G4double sval);
139  void ComputeLPMfunctions(G4double &fXiS, G4double &fGS, G4double &fPhiS,
140  const G4double eps, const G4double egamma,
141  const G4int izet);
142  struct LPMFuncs {
143  LPMFuncs() : fIsInitialized(false), fISDelta(100.), fSLimit(2.) {}
147  std::vector<G4double> fLPMFuncG;
148  std::vector<G4double> fLPMFuncPhi;
149  };
150 
151 
152 private:
153 
154  // hide assignment operator
155  G4PairProductionRelModel & operator=
156  (const G4PairProductionRelModel &right) = delete;
158 
159 protected:
160  static const G4int gMaxZet;
161  //
162  static const G4double gLPMconstant;
163  //
164  static const G4double gXGL[8];
165  static const G4double gWGL[8];
166  static const G4double gFelLowZet[8];
167  static const G4double gFinelLowZet[8];
168  //
169  static const G4double gXSecFactor;
171  //
172  static std::vector<ElementData*> gElementData;
174  //
177  //
179  //
185 };
186 
187 
188 
189 //
190 // Bethe screening functions for the elastic (coherent) scattering:
191 // Bethe's phi1, phi2 coherent screening functions were computed numerically
192 // by using (the universal) atomic form factors computed based on the Thomas-
193 // Fermi model of the atom (using numerical solution of the Thomas-Fermi
194 // screening function instead of Moliere's analytical approximation). The
195 // numerical results can be well approximated (better than Butcher & Messel
196 // especially near the delta=1 limit) by:
197 // ## if delta <= 1.4
198 // phi1(delta) = 20.806 - delta*(3.190 - 0.5710*delta)
199 // phi2(delta) = 20.234 - delta*(2.126 - 0.0903*delta)
200 // ## if delta > 1.4
201 // phi1(delta) = phi2(delta) = 21.0190 - 4.145*ln(delta + 0.958)
202 // with delta = 136mc^2kZ^{-1/3}/[E(Eg-E)] = 136Z^{-1/3}eps0/[eps(1-eps)] where
203 // Eg is the initial photon energy, E is the total energy transferred to one of
204 // the e-/e+ pair, eps0 = mc^2/Eg and eps = E/Eg.
205 
206 inline
208  G4double &phi2)
209 {
210  if (delta > 1.4) {
211  phi1 = 21.0190 - 4.145*G4Log(delta + 0.958);
212  phi2 = phi1;
213  } else {
214  phi1 = 20.806 - delta*(3.190 - 0.5710*delta);
215  phi2 = 20.234 - delta*(2.126 - 0.0903*delta);
216  }
217 }
218 
219 
220 // Compute the value of the screening function 3*PHI1(delta) - PHI2(delta):
222 {
223  return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
224  : 42.184 - delta*(7.444 - 1.623*delta);
225 }
226 
227 
228 // Compute the value of the screening function 1.5*PHI1(delta) +0.5*PHI2(delta):
230 {
231  return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
232  : 41.326 - delta*(5.848 - 0.902*delta);
233 }
234 
235 
236 // Same as ScreenFunction1 and ScreenFunction2 but computes them at once
238  G4double &f1, G4double &f2)
239 {
240  if (delta > 1.4) {
241  f1 = 42.038 - 8.29*G4Log(delta + 0.958);
242  f2 = f1;
243  } else {
244  f1 = 42.184 - delta*(7.444 - 1.623*delta);
245  f2 = 41.326 - delta*(5.848 - 0.902*delta);
246  }
247 }
248 
249 #endif
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ParticleDefinition * fThePositron
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
static std::vector< ElementData * > gElementData
static const G4double gFinelLowZet[8]
const char * p
Definition: xmltok.h:285
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ScreenFunction2(const G4double delta)
static const G4double gEgLPMActivation
void GetLPMFunctions(G4double &lpmGs, G4double &lpmPhis, const G4double sval)
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
static const G4double emax
Float_t f2
G4ParticleDefinition * fTheElectron
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4ParticleDefinition * fTheGamma
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static const G4double gLPMconstant
Float_t Z
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void ComputeLPMGsPhis(G4double &funcGS, G4double &funcPhiS, const G4double varShat)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
G4ParticleChangeForGamma * fParticleChange
void ComputePhi12(const G4double delta, G4double &phi1, G4double &phi2)
double A(double temperature)
static const G4double gXGL[8]
Float_t f1
static const G4double gWGL[8]
G4double ScreenFunction1(const G4double delta)
static const G4double gFelLowZet[8]
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
int G4int
Definition: G4Types.hh:78
Definition: G4Pow.hh:56
G4PairProductionRelModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheHeitlerLPM")
void ComputeLPMfunctions(G4double &fXiS, G4double &fGS, G4double &fPhiS, const G4double eps, const G4double egamma, const G4int izet)
static const G4double eps
#define DBL_MAX
Definition: templates.hh:83
G4double ComputeXSectionPerAtom(G4double gammaEnergy, G4double Z)