Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BetheHeitlerModel.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: G4BetheHeitlerModel
33//
34// Author: Vladimir Ivanchenko on base of Michel Maire code
35//
36// Creation date: 19.04.2005
37//
38// Modifications by Vladimir Ivanchenko, Michel Maire, Mihaly Novak
39//
40// Class Description:
41//
42// Implementation of gamma conversion to e+e- in the field of a nucleus
43// For details see Physics Reference Manual
44
45// -------------------------------------------------------------------
46//
47
48#ifndef G4BetheHeitlerModel_h
49#define G4BetheHeitlerModel_h 1
50
51#include "G4VEmModel.hh"
52#include "G4PhysicsTable.hh"
53#include "G4Log.hh"
54
55#include <vector>
56
58class G4Pow;
59class G4EmElementXS;
60
62{
63
64public:
65
66 explicit G4BetheHeitlerModel(const G4ParticleDefinition* p = nullptr,
67 const G4String& nam = "BetheHeitler");
68
69 ~G4BetheHeitlerModel() override;
70
71 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
72
74 G4VEmModel* masterModel) override;
75
77 G4double kinEnergy,
78 G4double Z,
79 G4double A=0.,
80 G4double cut=0.,
81 G4double emax=DBL_MAX) override;
82
83 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
85 const G4DynamicParticle*,
86 G4double tmin,
87 G4double maxEnergy) override;
88
89 // hide assignment operator
92
93protected:
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
102
107
108 static const G4int gMaxZet;
109
116
119
120 static std::vector<ElementData*> gElementData;
121};
122
123//
124// Bethe screening functions for the elastic (coherent) scattering:
125// Bethe's phi1, phi2 coherent screening functions were computed numerically
126// by using (the universal) atomic form factors computed based on the Thomas-
127// Fermi model of the atom (using numerical solution of the Thomas-Fermi
128// screening function instead of Moliere's analytical approximation). The
129// numerical results can be well approximated (better than Butcher & Messel
130// especially near the delta=1 limit) by:
131// ## if delta <= 1.4
132// phi1(delta) = 20.806 - delta*(3.190 - 0.5710*delta)
133// phi2(delta) = 20.234 - delta*(2.126 - 0.0903*delta)
134// ## if delta > 1.4
135// phi1(delta) = phi2(delta) = 21.0190 - 4.145*ln(delta + 0.958)
136// with delta = 136mc^2kZ^{-1/3}/[E(Eg-E)] = 136Z^{-1/3}eps0/[eps(1-eps)] where
137// Eg is the initial photon energy, E is the total energy transferred to one of
138// the e-/e+ pair, eps0 = mc^2/Eg and eps = E/Eg.
139
140// Compute the value of the screening function 3*PHI1(delta) - PHI2(delta):
142{
143 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
144 : 42.184 - delta*(7.444 - 1.623*delta);
145}
146
147// Compute the value of the screening function 1.5*PHI1(delta) +0.5*PHI2(delta):
149{
150 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
151 : 41.326 - delta*(5.848 - 0.902*delta);
152}
153
154// Same as ScreenFunction1 and ScreenFunction2 but computes them at once
156 G4double &f1, G4double &f2)
157{
158 if (delta > 1.4) {
159 f1 = 42.038 - 8.29*G4Log(delta + 0.958);
160 f2 = f1;
161 } else {
162 f1 = 42.184 - delta*(7.444 - 1.623*delta);
163 f2 = 41.326 - delta*(5.848 - 0.902*delta);
164 }
165}
166
167#endif
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4ParticleChangeForGamma * fParticleChange
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4BetheHeitlerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheHeitler")
const G4ParticleDefinition * fTheElectron
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
G4BetheHeitlerModel & operator=(const G4BetheHeitlerModel &right)=delete
G4double ScreenFunction1(const G4double delta)
static const G4int gMaxZet
G4BetheHeitlerModel(const G4BetheHeitlerModel &)=delete
static std::vector< ElementData * > gElementData
const G4ParticleDefinition * fTheGamma
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ParticleDefinition * fThePositron
G4double ScreenFunction2(const G4double delta)
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
Definition G4Pow.hh:49
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
#define DBL_MAX
Definition templates.hh:62