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.12 2008/10/15 15:54:57 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
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-20 SelectRandomElement is called for any initial gamma energy |
---|
46 | // in order to have selected element for polarized model (VI) |
---|
47 | // |
---|
48 | // Class Description: |
---|
49 | // |
---|
50 | // ------------------------------------------------------------------- |
---|
51 | // |
---|
52 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
53 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
54 | |
---|
55 | #include "G4BetheHeitlerModel.hh" |
---|
56 | #include "G4Electron.hh" |
---|
57 | #include "G4Positron.hh" |
---|
58 | #include "G4Gamma.hh" |
---|
59 | #include "Randomize.hh" |
---|
60 | #include "G4DataVector.hh" |
---|
61 | #include "G4PhysicsLogVector.hh" |
---|
62 | #include "G4ParticleChangeForGamma.hh" |
---|
63 | #include "G4LossTableManager.hh" |
---|
64 | |
---|
65 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
66 | |
---|
67 | using namespace std; |
---|
68 | |
---|
69 | G4BetheHeitlerModel::G4BetheHeitlerModel(const G4ParticleDefinition*, |
---|
70 | const G4String& nam) |
---|
71 | : G4VEmModel(nam), |
---|
72 | theCrossSectionTable(0), |
---|
73 | nbins(10) |
---|
74 | { |
---|
75 | fParticleChange = 0; |
---|
76 | theGamma = G4Gamma::Gamma(); |
---|
77 | thePositron = G4Positron::Positron(); |
---|
78 | theElectron = G4Electron::Electron(); |
---|
79 | } |
---|
80 | |
---|
81 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
82 | |
---|
83 | G4BetheHeitlerModel::~G4BetheHeitlerModel() |
---|
84 | { |
---|
85 | if(theCrossSectionTable) { |
---|
86 | theCrossSectionTable->clearAndDestroy(); |
---|
87 | delete theCrossSectionTable; |
---|
88 | } |
---|
89 | } |
---|
90 | |
---|
91 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
92 | |
---|
93 | void G4BetheHeitlerModel::Initialise(const G4ParticleDefinition*, |
---|
94 | const G4DataVector&) |
---|
95 | { |
---|
96 | if(!fParticleChange) { |
---|
97 | if(pParticleChange) { |
---|
98 | fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
99 | } else { |
---|
100 | fParticleChange = new G4ParticleChangeForGamma(); |
---|
101 | } |
---|
102 | } |
---|
103 | |
---|
104 | if(theCrossSectionTable) { |
---|
105 | theCrossSectionTable->clearAndDestroy(); |
---|
106 | delete theCrossSectionTable; |
---|
107 | } |
---|
108 | |
---|
109 | const G4ElementTable* theElementTable = G4Element::GetElementTable(); |
---|
110 | size_t nvect = G4Element::GetNumberOfElements(); |
---|
111 | theCrossSectionTable = new G4PhysicsTable(nvect); |
---|
112 | G4PhysicsLogVector* ptrVector; |
---|
113 | G4double emin = LowEnergyLimit(); |
---|
114 | G4double emax = HighEnergyLimit(); |
---|
115 | G4int n = nbins*G4int(log10(emax/emin)); |
---|
116 | G4bool spline = G4LossTableManager::Instance()->SplineFlag(); |
---|
117 | G4double e, value; |
---|
118 | |
---|
119 | for(size_t j=0; j<nvect ; j++) { |
---|
120 | |
---|
121 | ptrVector = new G4PhysicsLogVector(emin, emax, n); |
---|
122 | ptrVector->SetSpline(spline); |
---|
123 | G4double Z = (*theElementTable)[j]->GetZ(); |
---|
124 | G4int iz = G4int(Z); |
---|
125 | indexZ[iz] = j; |
---|
126 | |
---|
127 | for(G4int i=0; i<nbins; i++) { |
---|
128 | e = ptrVector->GetLowEdgeEnergy( i ) ; |
---|
129 | value = ComputeCrossSectionPerAtom(theGamma, e, Z); |
---|
130 | ptrVector->PutValue( i, value ); |
---|
131 | } |
---|
132 | |
---|
133 | theCrossSectionTable->insert(ptrVector); |
---|
134 | } |
---|
135 | } |
---|
136 | |
---|
137 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
138 | |
---|
139 | G4double G4BetheHeitlerModel::ComputeCrossSectionPerAtom( |
---|
140 | const G4ParticleDefinition*, |
---|
141 | G4double GammaEnergy, G4double Z, |
---|
142 | G4double, G4double, G4double) |
---|
143 | // Calculates the microscopic cross section in GEANT4 internal units. |
---|
144 | // A parametrized formula from L. Urban is used to estimate |
---|
145 | // the total cross section. |
---|
146 | // It gives a good description of the data from 1.5 MeV to 100 GeV. |
---|
147 | // below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass) |
---|
148 | // *(GammaEnergy-2electronmass) |
---|
149 | { |
---|
150 | static const G4double GammaEnergyLimit = 1.5*MeV; |
---|
151 | G4double CrossSection = 0.0 ; |
---|
152 | if ( Z < 1. ) return CrossSection; |
---|
153 | if ( GammaEnergy <= 2.0*electron_mass_c2 ) return CrossSection; |
---|
154 | |
---|
155 | static const G4double |
---|
156 | a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn, |
---|
157 | a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn; |
---|
158 | |
---|
159 | static const G4double |
---|
160 | b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381 *microbarn, |
---|
161 | b3= 1.3063 *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn; |
---|
162 | |
---|
163 | static const G4double |
---|
164 | c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn, |
---|
165 | c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn; |
---|
166 | |
---|
167 | G4double GammaEnergySave = GammaEnergy; |
---|
168 | if (GammaEnergy < GammaEnergyLimit) GammaEnergy = GammaEnergyLimit ; |
---|
169 | |
---|
170 | G4double X=log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X; |
---|
171 | |
---|
172 | G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5, |
---|
173 | F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5, |
---|
174 | F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5; |
---|
175 | |
---|
176 | CrossSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3); |
---|
177 | |
---|
178 | if (GammaEnergySave < GammaEnergyLimit) { |
---|
179 | |
---|
180 | X = (GammaEnergySave - 2.*electron_mass_c2) |
---|
181 | / (GammaEnergyLimit - 2.*electron_mass_c2); |
---|
182 | CrossSection *= X*X; |
---|
183 | } |
---|
184 | |
---|
185 | if (CrossSection < 0.) CrossSection = 0.; |
---|
186 | return CrossSection; |
---|
187 | } |
---|
188 | |
---|
189 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
190 | |
---|
191 | void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
192 | const G4MaterialCutsCouple* couple, |
---|
193 | const G4DynamicParticle* aDynamicGamma, |
---|
194 | G4double, |
---|
195 | G4double) |
---|
196 | // The secondaries e+e- energies are sampled using the Bethe - Heitler |
---|
197 | // cross sections with Coulomb correction. |
---|
198 | // A modified version of the random number techniques of Butcher & Messel |
---|
199 | // is used (Nuc Phys 20(1960),15). |
---|
200 | // |
---|
201 | // GEANT4 internal units. |
---|
202 | // |
---|
203 | // Note 1 : Effects due to the breakdown of the Born approximation at |
---|
204 | // low energy are ignored. |
---|
205 | // Note 2 : The differential cross section implicitly takes account of |
---|
206 | // pair creation in both nuclear and atomic electron fields. |
---|
207 | // However triplet prodution is not generated. |
---|
208 | { |
---|
209 | const G4Material* aMaterial = couple->GetMaterial(); |
---|
210 | |
---|
211 | G4double GammaEnergy = aDynamicGamma->GetKineticEnergy(); |
---|
212 | G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection(); |
---|
213 | |
---|
214 | G4double epsil ; |
---|
215 | G4double epsil0 = electron_mass_c2/GammaEnergy ; |
---|
216 | if(epsil0 > 1.0) return; |
---|
217 | |
---|
218 | // do it fast if GammaEnergy < 2. MeV |
---|
219 | static const G4double Egsmall=2.*MeV; |
---|
220 | |
---|
221 | // select randomly one element constituing the material |
---|
222 | const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy); |
---|
223 | |
---|
224 | if (GammaEnergy < Egsmall) { |
---|
225 | |
---|
226 | epsil = epsil0 + (0.5-epsil0)*G4UniformRand(); |
---|
227 | |
---|
228 | } else { |
---|
229 | // now comes the case with GammaEnergy >= 2. MeV |
---|
230 | |
---|
231 | // Extract Coulomb factor for this Element |
---|
232 | G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3()); |
---|
233 | if (GammaEnergy > 50.*MeV) FZ += 8.*(anElement->GetfCoulomb()); |
---|
234 | |
---|
235 | // limits of the screening variable |
---|
236 | G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3()); |
---|
237 | G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ; |
---|
238 | G4double screenmin = min(4.*screenfac,screenmax); |
---|
239 | |
---|
240 | // limits of the energy sampling |
---|
241 | G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ; |
---|
242 | G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin; |
---|
243 | |
---|
244 | // |
---|
245 | // sample the energy rate of the created electron (or positron) |
---|
246 | // |
---|
247 | //G4double epsil, screenvar, greject ; |
---|
248 | G4double screenvar, greject ; |
---|
249 | |
---|
250 | G4double F10 = ScreenFunction1(screenmin) - FZ; |
---|
251 | G4double F20 = ScreenFunction2(screenmin) - FZ; |
---|
252 | G4double NormF1 = max(F10*epsilrange*epsilrange,0.); |
---|
253 | G4double NormF2 = max(1.5*F20,0.); |
---|
254 | |
---|
255 | do { |
---|
256 | if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) { |
---|
257 | epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333); |
---|
258 | screenvar = screenfac/(epsil*(1-epsil)); |
---|
259 | greject = (ScreenFunction1(screenvar) - FZ)/F10; |
---|
260 | |
---|
261 | } else { |
---|
262 | epsil = epsilmin + epsilrange*G4UniformRand(); |
---|
263 | screenvar = screenfac/(epsil*(1-epsil)); |
---|
264 | greject = (ScreenFunction2(screenvar) - FZ)/F20; |
---|
265 | } |
---|
266 | |
---|
267 | } while( greject < G4UniformRand() ); |
---|
268 | |
---|
269 | } // end of epsil sampling |
---|
270 | |
---|
271 | // |
---|
272 | // fixe charges randomly |
---|
273 | // |
---|
274 | |
---|
275 | G4double ElectTotEnergy, PositTotEnergy; |
---|
276 | if (G4UniformRand() > 0.5) { |
---|
277 | |
---|
278 | ElectTotEnergy = (1.-epsil)*GammaEnergy; |
---|
279 | PositTotEnergy = epsil*GammaEnergy; |
---|
280 | |
---|
281 | } else { |
---|
282 | |
---|
283 | PositTotEnergy = (1.-epsil)*GammaEnergy; |
---|
284 | ElectTotEnergy = epsil*GammaEnergy; |
---|
285 | } |
---|
286 | |
---|
287 | // |
---|
288 | // scattered electron (positron) angles. ( Z - axis along the parent photon) |
---|
289 | // |
---|
290 | // universal distribution suggested by L. Urban |
---|
291 | // (Geant3 manual (1993) Phys211), |
---|
292 | // derived from Tsai distribution (Rev Mod Phys 49,421(1977)) |
---|
293 | |
---|
294 | G4double u; |
---|
295 | const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ; |
---|
296 | |
---|
297 | if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1; |
---|
298 | else u= - log(G4UniformRand()*G4UniformRand())/a2; |
---|
299 | |
---|
300 | G4double TetEl = u*electron_mass_c2/ElectTotEnergy; |
---|
301 | G4double TetPo = u*electron_mass_c2/PositTotEnergy; |
---|
302 | G4double Phi = twopi * G4UniformRand(); |
---|
303 | G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl); |
---|
304 | G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo); |
---|
305 | |
---|
306 | // |
---|
307 | // kinematic of the created pair |
---|
308 | // |
---|
309 | // the electron and positron are assumed to have a symetric |
---|
310 | // angular distribution with respect to the Z axis along the parent photon. |
---|
311 | |
---|
312 | G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2); |
---|
313 | |
---|
314 | G4ThreeVector ElectDirection (dxEl, dyEl, dzEl); |
---|
315 | ElectDirection.rotateUz(GammaDirection); |
---|
316 | |
---|
317 | // create G4DynamicParticle object for the particle1 |
---|
318 | G4DynamicParticle* aParticle1= new G4DynamicParticle( |
---|
319 | theElectron,ElectDirection,ElectKineEnergy); |
---|
320 | |
---|
321 | // the e+ is always created (even with Ekine=0) for further annihilation. |
---|
322 | |
---|
323 | G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2); |
---|
324 | |
---|
325 | G4ThreeVector PositDirection (dxPo, dyPo, dzPo); |
---|
326 | PositDirection.rotateUz(GammaDirection); |
---|
327 | |
---|
328 | // create G4DynamicParticle object for the particle2 |
---|
329 | G4DynamicParticle* aParticle2= new G4DynamicParticle( |
---|
330 | thePositron,PositDirection,PositKineEnergy); |
---|
331 | |
---|
332 | // Fill output vector |
---|
333 | fvect->push_back(aParticle1); |
---|
334 | fvect->push_back(aParticle2); |
---|
335 | |
---|
336 | // kill incident photon |
---|
337 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
338 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
339 | } |
---|
340 | |
---|
341 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
342 | |
---|
343 | |
---|