Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BetheHeitlerModel.cc
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 file
31//
32//
33// File name: G4BetheHeitlerModel
34//
35// Author: Vladimir Ivanchenko on base of Michel Maire code
36//
37// Creation date: 15.03.2005
38//
39// Modifications:
40// 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
41// 24-06-05 Increase number of bins to 200 (V.Ivantchenko)
42// 16-11-05 replace shootBit() by G4UniformRand() mma
43// 04-12-05 SetProposedKineticEnergy(0.) for the killed photon (mma)
44// 20-02-07 SelectRandomElement is called for any initial gamma energy
45// in order to have selected element for polarized model (VI)
46// 25-10-10 Removed unused table, added element selector (VI)
47//
48// Class Description:
49//
50// -------------------------------------------------------------------
51//
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
57#include "G4SystemOfUnits.hh"
58#include "G4Electron.hh"
59#include "G4Positron.hh"
60#include "G4Gamma.hh"
61#include "Randomize.hh"
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66using namespace std;
67
69 const G4String& nam)
70 : G4VEmModel(nam)
71{
72 fParticleChange = 0;
73 theGamma = G4Gamma::Gamma();
74 thePositron = G4Positron::Positron();
75 theElectron = G4Electron::Electron();
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81{}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86 const G4DataVector& cuts)
87{
88 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
96 G4double GammaEnergy, G4double Z,
98// Calculates the microscopic cross section in GEANT4 internal units.
99// A parametrized formula from L. Urban is used to estimate
100// the total cross section.
101// It gives a good description of the data from 1.5 MeV to 100 GeV.
102// below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
103// *(GammaEnergy-2electronmass)
104{
105 static const G4double GammaEnergyLimit = 1.5*MeV;
106 G4double xSection = 0.0 ;
107 if ( Z < 0.9 || GammaEnergy <= 2.0*electron_mass_c2 ) { return xSection; }
108
109 static const G4double
110 a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn,
111 a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn;
112
113 static const G4double
114 b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381 *microbarn,
115 b3= 1.3063 *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn;
116
117 static const G4double
118 c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn,
119 c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn;
120
121 G4double GammaEnergySave = GammaEnergy;
122 if (GammaEnergy < GammaEnergyLimit) { GammaEnergy = GammaEnergyLimit; }
123
124 G4double X=log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
125
126 G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
127 F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
128 F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;
129
130 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
131
132 if (GammaEnergySave < GammaEnergyLimit) {
133
134 X = (GammaEnergySave - 2.*electron_mass_c2)
135 / (GammaEnergyLimit - 2.*electron_mass_c2);
136 xSection *= X*X;
137 }
138
139 if (xSection < 0.) { xSection = 0.; }
140 return xSection;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144
145void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
146 const G4MaterialCutsCouple* couple,
147 const G4DynamicParticle* aDynamicGamma,
148 G4double,
149 G4double)
150// The secondaries e+e- energies are sampled using the Bethe - Heitler
151// cross sections with Coulomb correction.
152// A modified version of the random number techniques of Butcher & Messel
153// is used (Nuc Phys 20(1960),15).
154//
155// GEANT4 internal units.
156//
157// Note 1 : Effects due to the breakdown of the Born approximation at
158// low energy are ignored.
159// Note 2 : The differential cross section implicitly takes account of
160// pair creation in both nuclear and atomic electron fields.
161// However triplet prodution is not generated.
162{
163 const G4Material* aMaterial = couple->GetMaterial();
164
165 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
166 G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
167
168 G4double epsil ;
169 G4double epsil0 = electron_mass_c2/GammaEnergy ;
170 if(epsil0 > 1.0) { return; }
171
172 // do it fast if GammaEnergy < 2. MeV
173 static const G4double Egsmall=2.*MeV;
174
175 // select randomly one element constituing the material
176 const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
177
178 if (GammaEnergy < Egsmall) {
179
180 epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
181
182 } else {
183 // now comes the case with GammaEnergy >= 2. MeV
184
185 // Extract Coulomb factor for this Element
186 G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
187 if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
188
189 // limits of the screening variable
190 G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
191 G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
192 G4double screenmin = min(4.*screenfac,screenmax);
193
194 // limits of the energy sampling
195 G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
196 G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
197
198 //
199 // sample the energy rate of the created electron (or positron)
200 //
201 //G4double epsil, screenvar, greject ;
202 G4double screenvar, greject ;
203
204 G4double F10 = ScreenFunction1(screenmin) - FZ;
205 G4double F20 = ScreenFunction2(screenmin) - FZ;
206 G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
207 G4double NormF2 = max(1.5*F20,0.);
208
209 do {
210 if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
211 epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333);
212 screenvar = screenfac/(epsil*(1-epsil));
213 greject = (ScreenFunction1(screenvar) - FZ)/F10;
214
215 } else {
216 epsil = epsilmin + epsilrange*G4UniformRand();
217 screenvar = screenfac/(epsil*(1-epsil));
218 greject = (ScreenFunction2(screenvar) - FZ)/F20;
219 }
220
221 } while( greject < G4UniformRand() );
222
223 } // end of epsil sampling
224
225 //
226 // fixe charges randomly
227 //
228
229 G4double ElectTotEnergy, PositTotEnergy;
230 if (G4UniformRand() > 0.5) {
231
232 ElectTotEnergy = (1.-epsil)*GammaEnergy;
233 PositTotEnergy = epsil*GammaEnergy;
234
235 } else {
236
237 PositTotEnergy = (1.-epsil)*GammaEnergy;
238 ElectTotEnergy = epsil*GammaEnergy;
239 }
240
241 //
242 // scattered electron (positron) angles. ( Z - axis along the parent photon)
243 //
244 // universal distribution suggested by L. Urban
245 // (Geant3 manual (1993) Phys211),
246 // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
247
248 G4double u;
249 const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
250
251 if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
252 else u= - log(G4UniformRand()*G4UniformRand())/a2;
253
254 G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
255 G4double TetPo = u*electron_mass_c2/PositTotEnergy;
256 G4double Phi = twopi * G4UniformRand();
257 G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
258 G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
259
260 //
261 // kinematic of the created pair
262 //
263 // the electron and positron are assumed to have a symetric
264 // angular distribution with respect to the Z axis along the parent photon.
265
266 G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
267
268 G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
269 ElectDirection.rotateUz(GammaDirection);
270
271 // create G4DynamicParticle object for the particle1
272 G4DynamicParticle* aParticle1= new G4DynamicParticle(
273 theElectron,ElectDirection,ElectKineEnergy);
274
275 // the e+ is always created (even with Ekine=0) for further annihilation.
276
277 G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
278
279 G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
280 PositDirection.rotateUz(GammaDirection);
281
282 // create G4DynamicParticle object for the particle2
283 G4DynamicParticle* aParticle2= new G4DynamicParticle(
284 thePositron,PositDirection,PositKineEnergy);
285
286 // Fill output vector
287 fvect->push_back(aParticle1);
288 fvect->push_back(aParticle2);
289
290 // kill incident photon
291 fParticleChange->SetProposedKineticEnergy(0.);
292 fParticleChange->ProposeTrackStatus(fStopAndKill);
293}
294
295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
#define F10
#define F20
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
G4BetheHeitlerModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:209
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetlogZ3() const
G4double GetZ3() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)