Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BetheHeitlerModel.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: G4BetheHeitlerModel.hh 110527 2018-05-29 06:09:58Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4BetheHeitlerModel
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 19.04.2005
38 //
39 // Modifications:
40 // 02-02-06 Remove InitialiseCrossSectionPerAtom();
41 // 28-05-18 New version with improved screening function approximation, improved
42 // efficiency, documentation and cleanup (M. Novak)
43 //
44 // Class Description:
45 //
46 // Implementation of gamma convertion to e+e- in the field of a nucleus
47 //
48 
49 // -------------------------------------------------------------------
50 //
51 
52 #ifndef G4BetheHeitlerModel_h
53 #define G4BetheHeitlerModel_h 1
54 
55 #include "G4VEmModel.hh"
56 #include "G4PhysicsTable.hh"
57 #include "G4Log.hh"
58 
59 #include <vector>
60 
62 class G4Pow;
63 
65 {
66 
67 public:
68 
69  explicit G4BetheHeitlerModel(const G4ParticleDefinition* p = 0,
70  const G4String& nam = "BetheHeitler");
71 
72  virtual ~G4BetheHeitlerModel();
73 
74  virtual void Initialise(const G4ParticleDefinition*,
75  const G4DataVector&) override;
76 
77  virtual void InitialiseLocal(const G4ParticleDefinition*,
78  G4VEmModel* masterModel) override;
79 
81  G4double kinEnergy,
82  G4double Z,
83  G4double A=0.,
84  G4double cut=0.,
85  G4double emax=DBL_MAX) override;
86 
87  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
88  const G4MaterialCutsCouple*,
89  const G4DynamicParticle*,
90  G4double tmin,
91  G4double maxEnergy) override;
92 
93 protected:
94 
95  inline G4double ScreenFunction1(const G4double delta);
96 
97  inline G4double ScreenFunction2(const G4double delta);
98 
99  inline void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2);
100 
101  void InitialiseElementData();
102 
103  struct ElementData {
106  };
107 
108 private:
109 
110  // hide assignment operator
112  G4BetheHeitlerModel(const G4BetheHeitlerModel&) = delete;
113 
114 protected:
115 
116  static const G4int gMaxZet;
117 
123 
124  static std::vector<ElementData*> gElementData;
125 };
126 
127 
128 //
129 // Bethe screening functions for the elastic (coherent) scattering:
130 // Bethe's phi1, phi2 coherent screening functions were computed numerically
131 // by using (the universal) atomic form factors computed based on the Thomas-
132 // Fermi model of the atom (using numerical solution of the Thomas-Fermi
133 // screening function instead of Moliere's analytical approximation). The
134 // numerical results can be well approximated (better than Butcher & Messel
135 // especially near the delta=1 limit) by:
136 // ## if delta <= 1.4
137 // phi1(delta) = 20.806 - delta*(3.190 - 0.5710*delta)
138 // phi2(delta) = 20.234 - delta*(2.126 - 0.0903*delta)
139 // ## if delta > 1.4
140 // phi1(delta) = phi2(delta) = 21.0190 - 4.145*ln(delta + 0.958)
141 // with delta = 136mc^2kZ^{-1/3}/[E(Eg-E)] = 136Z^{-1/3}eps0/[eps(1-eps)] where
142 // Eg is the initial photon energy, E is the total energy transferred to one of
143 // the e-/e+ pair, eps0 = mc^2/Eg and eps = E/Eg.
144 
145 // Compute the value of the screening function 3*PHI1(delta) - PHI2(delta):
147 {
148  return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
149  : 42.184 - delta*(7.444 - 1.623*delta);
150 }
151 
152 
153 // Compute the value of the screening function 1.5*PHI1(delta) +0.5*PHI2(delta):
155 {
156  return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
157  : 41.326 - delta*(5.848 - 0.902*delta);
158 }
159 
160 
161 // Same as ScreenFunction1 and ScreenFunction2 but computes them at once
163  G4double &f1, G4double &f2)
164 {
165  if (delta > 1.4) {
166  f1 = 42.038 - 8.29*G4Log(delta + 0.958);
167  f2 = f1;
168  } else {
169  f1 = 42.184 - delta*(7.444 - 1.623*delta);
170  f2 = 41.326 - delta*(5.848 - 0.902*delta);
171  }
172 }
173 
174 
175 #endif
G4double ScreenFunction2(const G4double delta)
G4double ScreenFunction1(const G4double delta)
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
G4ParticleDefinition * fTheGamma
const char * p
Definition: xmltok.h:285
static std::vector< ElementData * > gElementData
G4BetheHeitlerModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
static const G4double emax
Float_t f2
static const G4int gMaxZet
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
double A(double temperature)
Float_t f1
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4BetheHeitlerModel & operator=(const G4BetheHeitlerModel &right)=delete
G4ParticleChangeForGamma * fParticleChange
int G4int
Definition: G4Types.hh:78
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
Definition: G4Pow.hh:56
G4ParticleDefinition * fThePositron
G4ParticleDefinition * fTheElectron
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
#define DBL_MAX
Definition: templates.hh:83