source: trunk/source/processes/electromagnetic/standard/src/G4BetheHeitlerModel.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

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