144{
145 const G4double sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
146 const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
147 Bohr_radius*Bohr_radius/(hbarc*hbarc);
148 const G4double epsmin = 1.e-4 , epsmax = 1.e10;
149
150 const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38., 47.,
151 50., 56., 64., 74., 79., 82. };
152
153 const G4double Tdat[22] = { 100*eV, 200*eV, 400*eV, 700*eV,
154 1*keV, 2*keV, 4*keV, 7*keV,
155 10*keV, 20*keV, 40*keV, 70*keV,
156 100*keV, 200*keV, 400*keV, 700*keV,
157 1*MeV, 2*MeV, 4*MeV, 7*MeV,
158 10*MeV, 20*MeV};
159
160
162 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
163 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
164 1.112,1.108,1.100,1.093,1.089,1.087 },
165 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
166 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
167 1.109,1.105,1.097,1.090,1.086,1.082 },
168 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
169 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
170 1.131,1.124,1.113,1.104,1.099,1.098 },
171 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
172 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
173 1.112,1.105,1.096,1.089,1.085,1.098 },
174 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
175 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
176 1.073,1.070,1.064,1.059,1.056,1.056 },
177 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
178 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
179 1.074,1.070,1.063,1.059,1.056,1.052 },
180 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
181 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
182 1.068,1.064,1.059,1.054,1.051,1.050 },
183 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
184 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
185 1.039,1.037,1.034,1.031,1.030,1.036 },
186 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
187 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
188 1.031,1.028,1.024,1.022,1.021,1.024 },
189 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
190 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
191 1.020,1.017,1.015,1.013,1.013,1.020 },
192 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
193 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
194 0.995,0.993,0.993,0.993,0.993,1.011 },
195 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
196 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
197 0.974,0.972,0.973,0.974,0.975,0.987 },
198 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
199 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
200 0.950,0.947,0.949,0.952,0.954,0.963 },
201 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
202 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
203 0.941,0.938,0.940,0.944,0.946,0.954 },
204 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
205 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
206 0.933,0.930,0.933,0.936,0.939,0.949 }};
207
209 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
210 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
211 1.131,1.126,1.117,1.108,1.103,1.100 },
212 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
213 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
214 1.138,1.132,1.122,1.113,1.108,1.102 },
215 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
216 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
217 1.203,1.190,1.173,1.159,1.151,1.145 },
218 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
219 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
220 1.225,1.210,1.191,1.175,1.166,1.174 },
221 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
222 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
223 1.217,1.203,1.184,1.169,1.160,1.151 },
224 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
225 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
226 1.237,1.222,1.201,1.184,1.174,1.159 },
227 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
228 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
229 1.252,1.234,1.212,1.194,1.183,1.170 },
230 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
231 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
232 1.254,1.237,1.214,1.195,1.185,1.179 },
233 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
234 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
235 1.312,1.288,1.258,1.235,1.221,1.205 },
236 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
237 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
238 1.320,1.294,1.264,1.240,1.226,1.214 },
239 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
240 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
241 1.328,1.302,1.270,1.245,1.231,1.233 },
242 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
243 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
244 1.361,1.330,1.294,1.267,1.251,1.239 },
245 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
246 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
247 1.409,1.372,1.330,1.298,1.280,1.258 },
248 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
249 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
250 1.442,1.400,1.354,1.319,1.299,1.272 },
251 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
252 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
253 1.456,1.412,1.364,1.328,1.307,1.282 }};
254
255
257 G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
258 ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
259 G4double bg2lim = Tlim*(Tlim+2.*electron_mass_c2)/
260 (electron_mass_c2*electron_mass_c2);
261
262 G4double sig0[15] = {0.2672*barn, 0.5922*barn, 2.653*barn, 6.235*barn,
263 11.69*barn , 13.24*barn , 16.12*barn, 23.00*barn ,
264 35.13*barn , 39.95*barn , 50.85*barn, 67.19*barn ,
265 91.15*barn , 104.4*barn , 113.1*barn};
266
267 G4double hecorr[15] = {120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
268 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
269 -22.30};
270
272 SetParticle(part);
273
274 G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
275
276
277
278
279
280 G4double eKineticEnergy = KineticEnergy;
281
282 if(mass > electron_mass_c2)
283 {
285 G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
287 G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
288 eKineticEnergy = electron_mass_c2*tau ;
289 }
290
291 G4double ChargeSquare = charge*charge;
292
293 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
294 G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
295 /(eTotalEnergy*eTotalEnergy);
296 G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
297 /(electron_mass_c2*electron_mass_c2);
298
300
301 if (eps<epsmin) sigma = 2.*eps*eps;
302 else if(eps<epsmax) sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
303 else sigma = log(2.*eps)-1.+1./eps;
304
305 sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
306
307
309
310
312 while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
313 if (iZ==14) iZ = 13;
314 if (iZ==-1) iZ = 0 ;
315
318 G4double ratZ = (AtomicNumber-Z1)*(AtomicNumber+Z1)/
319 ((Z2-Z1)*(Z2+Z1));
320
321 if(eKineticEnergy <= Tlim)
322 {
323
325 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
326 if(iT==21) iT = 20;
327 if(iT==-1) iT = 0 ;
328
329
330 G4double T = Tdat[iT], E = T + electron_mass_c2;
331 G4double b2small = T*(E+electron_mass_c2)/(E*E);
332
333 T = Tdat[iT+1]; E = T + electron_mass_c2;
334 G4double b2big = T*(E+electron_mass_c2)/(E*E);
335 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
336
337 if (charge < 0.)
338 {
339 c1 = celectron[iZ][iT];
340 c2 = celectron[iZ+1][iT];
341 cc1 = c1+ratZ*(c2-c1);
342
343 c1 = celectron[iZ][iT+1];
344 c2 = celectron[iZ+1][iT+1];
345 cc2 = c1+ratZ*(c2-c1);
346
347 corr = cc1+ratb2*(cc2-cc1);
348
349 sigma *= sigmafactor/corr;
350 }
351 else
352 {
353 c1 = cpositron[iZ][iT];
354 c2 = cpositron[iZ+1][iT];
355 cc1 = c1+ratZ*(c2-c1);
356
357 c1 = cpositron[iZ][iT+1];
358 c2 = cpositron[iZ+1][iT+1];
359 cc2 = c1+ratZ*(c2-c1);
360
361 corr = cc1+ratb2*(cc2-cc1);
362
363 sigma *= sigmafactor/corr;
364 }
365 }
366 else
367 {
368 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
369 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
370 if((AtomicNumber >= Z1) && (AtomicNumber <= Z2))
371 sigma = c1+ratZ*(c2-c1) ;
372 else if(AtomicNumber < Z1)
373 sigma = AtomicNumber*AtomicNumber*c1/(Z1*Z1);
374 else if(AtomicNumber > Z2)
375 sigma = AtomicNumber*AtomicNumber*c2/(Z2*Z2);
376 }
377 return sigma;
378
379}