Geant4 11.4.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4PairProductionRelModel
33//
34// Author: Andreas Schaelicke
35//
36// Creation date: 02.04.2009
37//
38// Modifications:
39// 28-05-18 New version with improved screening function approximation, improved
40// LPM function approximation, efficiency, documentation and cleanup.
41// Corrected call to selecting target atom in the final state sampling.
42// (M. Novak)
43//
44// Class Description:
45//
46// Implementation of gamma conversion to e+e- in the field of a nucleus
47// relativistic approximation
48//
49
50// -------------------------------------------------------------------
51//
52
53#ifndef G4PairProductionRelModel_h
54#define G4PairProductionRelModel_h 1
55
57
58#include "G4VEmModel.hh"
59#include "G4Log.hh"
60#include <vector>
61
63class G4Pow;
64
66{
67
68public:
69
70 explicit G4PairProductionRelModel(const G4ParticleDefinition* p = nullptr,
71 const G4String& nam = "BetheHeitlerLPM");
72
74
75 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
76
78 G4VEmModel* masterModel) override;
79
81 G4double kinEnergy,
83 G4double A=0.,
84 G4double cut=0.,
85 G4double emax=DBL_MAX) override;
87 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89 const G4DynamicParticle*,
90 G4double tmin,
91 G4double maxEnergy) override;
92
94 const G4Material*,G4double) override;
95
96 inline void SetLPMflag(G4bool val) { fIsUseLPMCorrection = val; }
97 inline G4bool LPMflag() const { return fIsUseLPMCorrection; }
98
100 (const G4PairProductionRelModel &right) = delete;
102
103protected:
104
105 // for evaluating screening related functions
106 inline void ComputePhi12(const G4double delta,
107 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,
111 G4double &f1, G4double &f2);
112 // helper methods for cross-section computation under different approximations
115 G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy,
116 G4double Z);
118 G4double gammaEnergy, G4double Z);
119
120private:
121
122 // for creating some data structure per Z
123 void InitialiseElementData();
124 struct ElementData {
125 G4double fLogZ13;
126 G4double fCoulomb;
127 G4double fLradEl;
128 G4double fDeltaFactor;
129 G4double fDeltaMaxLow;
130 G4double fDeltaMaxHigh;
131 G4double fEtaValue;
132 G4double fLPMVarS1Cond;
133 G4double fLPMILVarS1Cond;
134 };
135 void ComputeLPMfunctions(G4double &fXiS, G4double &fGS, G4double &fPhiS,
136 const G4double eps, const G4double egamma,
137 const G4int izet);
138
139protected:
140 static const G4int gMaxZet;
141 //
142 static const G4double gLPMconstant;
143 //
144 static const G4double gXGL[8];
145 static const G4double gWGL[8];
146 static const G4double gFelLowZet[8];
147 static const G4double gFinelLowZet[8];
148 //
149 static const G4double gXSecFactor;
151 //
152 static std::vector<ElementData*> gElementData;
153 //
157 //
159 //
162 //
168};
169//
170// Bethe screening functions for the elastic (coherent) scattering:
171// Bethe's phi1, phi2 coherent screening functions were computed numerically
172// by using (the universal) atomic form factors computed based on the Thomas-
173// Fermi model of the atom (using numerical solution of the Thomas-Fermi
174// screening function instead of Moliere's analytical approximation). The
175// numerical results can be well approximated (better than Butcher & Messel
176// especially near the delta=1 limit) by:
177// ## if delta <= 1.4
178// phi1(delta) = 20.806 - delta*(3.190 - 0.5710*delta)
179// phi2(delta) = 20.234 - delta*(2.126 - 0.0903*delta)
180// ## if delta > 1.4
181// phi1(delta) = phi2(delta) = 21.0190 - 4.145*ln(delta + 0.958)
182// with delta = 136mc^2kZ^{-1/3}/[E(Eg-E)] = 136Z^{-1/3}eps0/[eps(1-eps)] where
183// Eg is the initial photon energy, E is the total energy transferred to one of
184// the e-/e+ pair, eps0 = mc^2/Eg and eps = E/Eg.
185
187 G4double &phi1,
188 G4double &phi2)
189{
190 if (delta > 1.4) {
191 phi1 = 21.0190 - 4.145*G4Log(delta + 0.958);
192 phi2 = phi1;
193 } else {
194 phi1 = 20.806 - delta*(3.190 - 0.5710*delta);
195 phi2 = 20.234 - delta*(2.126 - 0.0903*delta);
196 }
197}
198
199// Compute the value of the screening function 3*PHI1(delta) - PHI2(delta):
201{
202 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
203 : 42.184 - delta*(7.444 - 1.623*delta);
204}
205
206// Compute the value of the screening function 1.5*PHI1(delta) +0.5*PHI2(delta):
208{
209 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
210 : 41.326 - delta*(5.848 - 0.902*delta);
211}
212
213// Same as ScreenFunction1 and ScreenFunction2 but computes them at once
215 G4double &f1, G4double &f2)
216{
217 if (delta > 1.4) {
218 f1 = 42.038 - 8.29*G4Log(delta + 0.958);
219 f2 = f1;
220 } else {
221 f1 = 42.184 - delta*(7.444 - 1.623*delta);
222 f2 = 41.326 - delta*(5.848 - 0.902*delta);
223 }
224}
225
226#endif
G4double G4Log(G4double x)
Definition G4Log.hh:169
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
static const G4double gFelLowZet[8]
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeParametrizedXSectionPerAtom(G4double gammaEnergy, G4double Z)
static const G4double gEgLPMActivation
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
static const G4double gFinelLowZet[8]
static std::vector< ElementData * > gElementData
G4PairProductionRelModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheHeitlerLPM")
static const G4double gXGL[8]
G4double ScreenFunction1(const G4double delta)
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
G4ParticleDefinition * fThePositron
void ComputePhi12(const G4double delta, G4double &phi1, G4double &phi2)
G4double ComputeXSectionPerAtom(G4double gammaEnergy, G4double Z)
G4PairProductionRelModel(const G4PairProductionRelModel &)=delete
G4double ScreenFunction2(const G4double delta)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static const G4double gWGL[8]
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * fTheElectron
Definition G4Pow.hh:49
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
#define DBL_MAX
Definition templates.hh:62