source: trunk/source/processes/electromagnetic/standard/src/G4eBremsstrahlungModel.cc@ 1005

Last change on this file since 1005 was 991, checked in by garnier, 17 years ago

update

File size: 31.1 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: G4eBremsstrahlungModel.cc,v 1.43 2008/11/13 19:28:58 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4eBremsstrahlungModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 11-11-02 Fix division by 0 (V.Ivanchenko)
43// 04-12-02 Change G4DynamicParticle constructor in PostStep (V.Ivanchenko)
44// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
45// 24-01-03 Fix for compounds (V.Ivanchenko)
46// 27-01-03 Make models region aware (V.Ivanchenko)
47// 13-02-03 Add name (V.Ivanchenko)
48// 09-05-03 Fix problem of supression function + optimise sampling (V.Ivanchenko)
49// 20-05-04 Correction to ensure unit independence (L.Urban)
50// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
51// 03-08-05 Add extra protection at initialisation (V.Ivantchenko)
52// 07-02-06 public function ComputeCrossSectionPerAtom() (mma)
53// 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI)
54// 27-03-06 Fix calculation of fl parameter at low energy (energy loss) (VI)
55// 15-02-07 correct LPMconstant by a factor 2, thanks to G. Depaola (mma)
56// 09-09-08 MigdalConstant increased in (2pi)^2 times (A.Schaelicke)
57//
58// Class Description:
59//
60//
61// -------------------------------------------------------------------
62//
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66#include "G4eBremsstrahlungModel.hh"
67#include "G4Electron.hh"
68#include "G4Positron.hh"
69#include "G4Gamma.hh"
70#include "Randomize.hh"
71#include "G4Material.hh"
72#include "G4Element.hh"
73#include "G4ElementVector.hh"
74#include "G4ProductionCutsTable.hh"
75#include "G4DataVector.hh"
76#include "G4ParticleChangeForLoss.hh"
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
80using namespace std;
81
82G4eBremsstrahlungModel::G4eBremsstrahlungModel(const G4ParticleDefinition* p,
83 const G4String& nam)
84 : G4VEmModel(nam),
85 particle(0),
86 isElectron(true),
87 probsup(1.0),
88 MigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
89 LPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)),
90 isInitialised(false)
91{
92 if(p) SetParticle(p);
93 theGamma = G4Gamma::Gamma();
94 minThreshold = 1.0*keV;
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
99G4eBremsstrahlungModel::~G4eBremsstrahlungModel()
100{
101 size_t n = partialSumSigma.size();
102 if(n > 0) {
103 for(size_t i=0; i<n; i++) {
104 delete partialSumSigma[i];
105 }
106 }
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110
111void G4eBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p)
112{
113 particle = p;
114 if(p == G4Electron::Electron()) isElectron = true;
115 else isElectron = false;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120G4double G4eBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*,
121 const G4MaterialCutsCouple*)
122{
123 return minThreshold;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
128void G4eBremsstrahlungModel::Initialise(const G4ParticleDefinition* p,
129 const G4DataVector& cuts)
130{
131 if(p) SetParticle(p);
132 highKinEnergy = HighEnergyLimit();
133 lowKinEnergy = LowEnergyLimit();
134 const G4ProductionCutsTable* theCoupleTable=
135 G4ProductionCutsTable::GetProductionCutsTable();
136
137 if(theCoupleTable) {
138 G4int numOfCouples = theCoupleTable->GetTableSize();
139
140 G4int nn = partialSumSigma.size();
141 G4int nc = cuts.size();
142 if(nn > 0) {
143 for (G4int ii=0; ii<nn; ii++){
144 G4DataVector* a=partialSumSigma[ii];
145 if ( a ) delete a;
146 }
147 partialSumSigma.clear();
148 }
149 if(numOfCouples>0) {
150 for (G4int i=0; i<numOfCouples; i++) {
151 G4double cute = DBL_MAX;
152 if(i < nc) cute = cuts[i];
153 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
154 const G4Material* material = couple->GetMaterial();
155 G4DataVector* dv = ComputePartialSumSigma(material, 0.5*highKinEnergy,
156 std::min(cute, 0.25*highKinEnergy));
157 partialSumSigma.push_back(dv);
158 }
159 }
160 }
161 if(isInitialised) return;
162
163 if(pParticleChange) {
164 fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
165 } else {
166 fParticleChange = new G4ParticleChangeForLoss();
167 }
168 isInitialised = true;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
173G4double G4eBremsstrahlungModel::ComputeDEDXPerVolume(
174 const G4Material* material,
175 const G4ParticleDefinition* p,
176 G4double kineticEnergy,
177 G4double cutEnergy)
178{
179 if(!particle) SetParticle(p);
180 if(kineticEnergy < lowKinEnergy) return 0.0;
181
182 const G4double thigh = 100.*GeV;
183
184 G4double cut = std::min(cutEnergy, kineticEnergy);
185
186 G4double rate, loss;
187 const G4double factorHigh = 36./(1450.*GeV);
188 const G4double coef1 = -0.5;
189 const G4double coef2 = 2./9.;
190
191 const G4ElementVector* theElementVector = material->GetElementVector();
192 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
193
194 G4double totalEnergy = kineticEnergy + electron_mass_c2;
195 G4double dedx = 0.0;
196
197 // loop for elements in the material
198 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
199
200 G4double Z = (*theElementVector)[i]->GetZ();
201 G4double natom = theAtomicNumDensityVector[i];
202
203 // loss for MinKinEnergy<KineticEnergy<=100 GeV
204 if (kineticEnergy <= thigh) {
205
206 // x = log(totalEnergy/electron_mass_c2);
207 loss = ComputeBremLoss(Z, kineticEnergy, cut) ;
208 if (!isElectron) loss *= PositronCorrFactorLoss(Z, kineticEnergy, cut);
209
210 // extrapolation for KineticEnergy>100 GeV
211 } else {
212
213 // G4double xhigh = log(thigh/electron_mass_c2);
214 G4double cuthigh = thigh*0.5;
215
216 if (cut < thigh) {
217
218 loss = ComputeBremLoss(Z, thigh, cut) ;
219 if (!isElectron) loss *= PositronCorrFactorLoss(Z, thigh, cut) ;
220 rate = cut/totalEnergy;
221 loss *= (1. + coef1*rate + coef2*rate*rate);
222 rate = cut/thigh;
223 loss /= (1.+coef1*rate+coef2*rate*rate);
224
225 } else {
226
227 loss = ComputeBremLoss(Z, thigh, cuthigh) ;
228 if (!isElectron) loss *= PositronCorrFactorLoss(Z, thigh, cuthigh) ;
229 rate = cut/totalEnergy;
230 loss *= (1. + coef1*rate + coef2*rate*rate);
231 loss *= cut*factorHigh;
232 }
233 }
234 loss *= natom;
235
236 G4double kp2 = MigdalConstant*totalEnergy*totalEnergy
237 * (material->GetElectronDensity()) ;
238
239 // now compute the correction due to the supression(s)
240 G4double kmin = 1.*eV;
241 G4double kmax = cut;
242
243 if (kmax > kmin) {
244
245 G4double floss = 0.;
246 G4int nmax = 100;
247
248 G4double vmin=log(kmin);
249 G4double vmax=log(kmax) ;
250 G4int nn = (G4int)(nmax*(vmax-vmin)/(log(highKinEnergy)-vmin)) ;
251 G4double u,fac,c,v,dv ;
252 if(nn > 0) {
253
254 dv = (vmax-vmin)/nn ;
255 v = vmin-dv ;
256
257 for(G4int n=0; n<=nn; n++) {
258
259 v += dv;
260 u = exp(v);
261 fac = u*SupressionFunction(material,kineticEnergy,u);
262 fac *= probsup*(u*u/(u*u+kp2))+1.-probsup;
263 if ((n==0)||(n==nn)) c=0.5;
264 else c=1. ;
265 fac *= c ;
266 floss += fac ;
267 }
268 floss *=dv/(kmax-kmin);
269
270 } else {
271 floss = 1.;
272 }
273 if(floss > 1.) floss = 1.;
274 // correct the loss
275 loss *= floss;
276 }
277 dedx += loss;
278 }
279 if(dedx < 0.) dedx = 0.;
280 return dedx;
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
284
285G4double G4eBremsstrahlungModel::ComputeBremLoss(G4double Z, G4double T,
286 G4double Cut)
287
288// compute loss due to soft brems
289{
290 static const G4double beta=1.0, ksi=2.0;
291 static const G4double clossh = 0.254 , closslow = 1./3. , alosslow = 1. ;
292 static const G4double Tlim= 10.*MeV ;
293
294 static const G4double xlim = 1.2 ;
295 static const G4int NZ = 8 ;
296 static const G4int Nloss = 11 ;
297 static const G4double ZZ[NZ] =
298 {2.,4.,6.,14.,26.,50.,82.,92.};
299 static const G4double coefloss[NZ][Nloss] = {
300 // Z=2
301 { 0.98916, 0.47564, -0.2505, -0.45186, 0.14462,
302 0.21307, -0.013738, -0.045689, -0.0042914, 0.0034429,
303 0.00064189},
304
305 // Z=4
306 { 1.0626, 0.37662, -0.23646, -0.45188, 0.14295,
307 0.22906, -0.011041, -0.051398, -0.0055123, 0.0039919,
308 0.00078003},
309 // Z=6
310 { 1.0954, 0.315, -0.24011, -0.43849, 0.15017,
311 0.23001, -0.012846, -0.052555, -0.0055114, 0.0041283,
312 0.00080318},
313
314 // Z=14
315 { 1.1649, 0.18976, -0.24972, -0.30124, 0.1555,
316 0.13565, -0.024765, -0.027047, -0.00059821, 0.0019373,
317 0.00027647},
318
319 // Z=26
320 { 1.2261, 0.14272, -0.25672, -0.28407, 0.13874,
321 0.13586, -0.020562, -0.026722, -0.00089557, 0.0018665,
322 0.00026981},
323
324 // Z=50
325 { 1.3147, 0.020049, -0.35543, -0.13927, 0.17666,
326 0.073746, -0.036076, -0.013407, 0.0025727, 0.00084005,
327 -1.4082e-05},
328
329 // Z=82
330 { 1.3986, -0.10586, -0.49187, -0.0048846, 0.23621,
331 0.031652, -0.052938, -0.0076639, 0.0048181, 0.00056486,
332 -0.00011995},
333
334 // Z=92
335 { 1.4217, -0.116, -0.55497, -0.044075, 0.27506,
336 0.081364, -0.058143, -0.023402, 0.0031322, 0.0020201,
337 0.00017519}
338
339 } ;
340 static G4double aaa = 0.414;
341 static G4double bbb = 0.345;
342 static G4double ccc = 0.460;
343
344 G4int iz = 0;
345 G4double delz = 1.e6;
346 for (G4int ii=0; ii<NZ; ii++)
347 {
348 G4double dz = std::abs(Z-ZZ[ii]);
349 if(dz < delz) {
350 iz = ii;
351 delz = dz;
352 }
353 }
354
355 G4double xx = log10(T/MeV);
356 G4double fl = 1.;
357
358 if (xx <= xlim)
359 {
360 xx /= xlim;
361 G4double yy = 1.0;
362 fl = 0.0;
363 for (G4int j=0; j<Nloss; j++) {
364 fl += yy+coefloss[iz][j];
365 yy *= xx;
366 }
367 if (fl < 0.00001) fl = 0.00001;
368 else if (fl > 1.0) fl = 1.0;
369 }
370
371 G4double loss;
372 G4double E = T+electron_mass_c2 ;
373
374 loss = Z*(Z+ksi)*E*E/(T+E)*exp(beta*log(Cut/T))*(2.-clossh*exp(log(Z)/4.));
375 if (T <= Tlim) loss /= exp(closslow*log(Tlim/T));
376 if( T <= Cut) loss *= exp(alosslow*log(T/Cut));
377 // correction
378 loss *= (aaa+bbb*T/Tlim)/(1.+ccc*T/Tlim);
379 loss *= fl;
380 loss /= Avogadro;
381
382 return loss;
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386
387G4double G4eBremsstrahlungModel::PositronCorrFactorLoss(G4double Z,
388 G4double kineticEnergy, G4double cut)
389
390//calculates the correction factor for the energy loss due to bremsstrahlung for positrons
391//the same correction is in the (discrete) bremsstrahlung
392
393{
394 static const G4double K = 132.9416*eV ;
395 static const G4double a1=4.15e-1, a3=2.10e-3, a5=54.0e-5 ;
396
397 G4double x = log(kineticEnergy/(K*Z*Z)), x2 = x*x, x3 = x2*x;
398 G4double eta = 0.5+atan(a1*x+a3*x3+a5*x3*x2)/pi;
399 G4double e0 = cut/kineticEnergy;
400
401 G4double factor = 0.0;
402 if (e0 < 1.0) {
403 factor=log(1.-e0)/eta;
404 factor=exp(factor);
405 }
406 factor = eta*(1.-factor)/e0;
407
408 return factor;
409}
410
411//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412
413G4double G4eBremsstrahlungModel::CrossSectionPerVolume(
414 const G4Material* material,
415 const G4ParticleDefinition* p,
416 G4double kineticEnergy,
417 G4double cutEnergy,
418 G4double maxEnergy)
419{
420 if(!particle) SetParticle(p);
421 G4double cross = 0.0;
422 G4double tmax = min(maxEnergy, kineticEnergy);
423 G4double cut = max(cutEnergy, minThreshold);
424 if(cut >= tmax) return cross;
425
426 const G4ElementVector* theElementVector = material->GetElementVector();
427 const G4double* theAtomNumDensityVector = material->GetAtomicNumDensityVector();
428 G4double dum=0.;
429
430 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
431
432 cross += theAtomNumDensityVector[i] * ComputeCrossSectionPerAtom(p,
433 kineticEnergy, (*theElementVector)[i]->GetZ(), dum, cut);
434 if (tmax < kineticEnergy) {
435 cross -= theAtomNumDensityVector[i] * ComputeCrossSectionPerAtom(p,
436 kineticEnergy, (*theElementVector)[i]->GetZ(), dum, tmax);
437 }
438 }
439
440 // now compute the correction due to the supression(s)
441
442 G4double kmax = tmax;
443 G4double kmin = cut;
444
445 G4double totalEnergy = kineticEnergy+electron_mass_c2 ;
446 G4double kp2 = MigdalConstant*totalEnergy*totalEnergy
447 *(material->GetElectronDensity());
448
449 G4double fsig = 0.;
450 G4int nmax = 100;
451 G4double vmin=log(kmin);
452 G4double vmax=log(kmax) ;
453 G4int nn = (G4int)(nmax*(vmax-vmin)/(log(highKinEnergy)-vmin));
454 G4double u,fac,c,v,dv,y ;
455 if(nn > 0) {
456
457 dv = (vmax-vmin)/nn ;
458 v = vmin-dv ;
459 for(G4int n=0; n<=nn; n++) {
460
461 v += dv;
462 u = exp(v);
463 fac = SupressionFunction(material, kineticEnergy, u);
464 y = u/kmax;
465 fac *= (4.-4.*y+3.*y*y)/3.;
466 fac *= probsup*(u*u/(u*u+kp2))+1.-probsup;
467
468 if ((n==0)||(n==nn)) c=0.5;
469 else c=1. ;
470
471 fac *= c;
472 fsig += fac;
473 }
474 y = kmin/kmax ;
475 fsig *=dv/(-4.*log(y)/3.-4.*(1.-y)/3.+0.5*(1.-y*y));
476
477 } else {
478
479 fsig = 1.;
480 }
481 if (fsig > 1.) fsig = 1.;
482
483 // correct the cross section
484 cross *= fsig;
485
486 return cross;
487}
488
489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
490
491G4double G4eBremsstrahlungModel::ComputeCrossSectionPerAtom(
492 const G4ParticleDefinition*,
493 G4double kineticEnergy,
494 G4double Z, G4double,
495 G4double cut, G4double)
496
497// Calculates the cross section per atom in GEANT4 internal units.
498//
499
500{
501 G4double cross = 0.0 ;
502 if ( kineticEnergy < 1*keV || kineticEnergy < cut) return cross;
503
504 static const G4double ksi=2.0, alfa=1.00;
505 static const G4double csigh = 0.127, csiglow = 0.25, asiglow = 0.020*MeV ;
506 static const G4double Tlim = 10.*MeV ;
507
508 static const G4double xlim = 1.2 ;
509 static const G4int NZ = 8 ;
510 static const G4int Nsig = 11 ;
511 static const G4double ZZ[NZ] =
512 {2.,4.,6.,14.,26.,50.,82.,92.} ;
513 static const G4double coefsig[NZ][Nsig] = {
514 // Z=2
515 { 0.4638, 0.37748, 0.32249, -0.060362, -0.065004,
516 -0.033457, -0.004583, 0.011954, 0.0030404, -0.0010077,
517 -0.00028131},
518
519 // Z=4
520 { 0.50008, 0.33483, 0.34364, -0.086262, -0.055361,
521 -0.028168, -0.0056172, 0.011129, 0.0027528, -0.00092265,
522 -0.00024348},
523
524 // Z=6
525 { 0.51587, 0.31095, 0.34996, -0.11623, -0.056167,
526 -0.0087154, 0.00053943, 0.0054092, 0.00077685, -0.00039635,
527 -6.7818e-05},
528
529 // Z=14
530 { 0.55058, 0.25629, 0.35854, -0.080656, -0.054308,
531 -0.049933, -0.00064246, 0.016597, 0.0021789, -0.001327,
532 -0.00025983},
533
534 // Z=26
535 { 0.5791, 0.26152, 0.38953, -0.17104, -0.099172,
536 0.024596, 0.023718, -0.0039205, -0.0036658, 0.00041749,
537 0.00023408},
538
539 // Z=50
540 { 0.62085, 0.27045, 0.39073, -0.37916, -0.18878,
541 0.23905, 0.095028, -0.068744, -0.023809, 0.0062408,
542 0.0020407},
543
544 // Z=82
545 { 0.66053, 0.24513, 0.35404, -0.47275, -0.22837,
546 0.35647, 0.13203, -0.1049, -0.034851, 0.0095046,
547 0.0030535},
548
549 // Z=92
550 { 0.67143, 0.23079, 0.32256, -0.46248, -0.20013,
551 0.3506, 0.11779, -0.1024, -0.032013, 0.0092279,
552 0.0028592}
553
554 } ;
555
556 G4int iz = 0 ;
557 G4double delz = 1.e6 ;
558 for (G4int ii=0; ii<NZ; ii++)
559 {
560 G4double absdelz = std::abs(Z-ZZ[ii]);
561 if(absdelz < delz)
562 {
563 iz = ii ;
564 delz = absdelz;
565 }
566 }
567
568 G4double xx = log10(kineticEnergy/MeV) ;
569 G4double fs = 1. ;
570
571 if (xx <= xlim) {
572
573 fs = coefsig[iz][Nsig-1] ;
574 for (G4int j=Nsig-2; j>=0; j--) {
575
576 fs = fs*xx+coefsig[iz][j] ;
577 }
578 if(fs < 0.) fs = 0.;
579 }
580
581 cross = Z*(Z+ksi)*(1.-csigh*exp(log(Z)/4.))*pow(log(kineticEnergy/cut),alfa);
582
583 if (kineticEnergy <= Tlim)
584 cross *= exp(csiglow*log(Tlim/kineticEnergy))
585 *(1.+asiglow/(sqrt(Z)*kineticEnergy));
586
587 if (!isElectron)
588 cross *= PositronCorrFactorSigma(Z, kineticEnergy, cut);
589
590 cross *= fs/Avogadro ;
591
592 if (cross < 0.) cross = 0.;
593 return cross;
594}
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
597
598G4double G4eBremsstrahlungModel::PositronCorrFactorSigma( G4double Z,
599 G4double kineticEnergy, G4double cut)
600
601//Calculates the correction factor for the total cross section of the positron
602// bremsstrahl.
603// Eta is the ratio of positron to electron energy loss by bremstrahlung.
604// A parametrized formula from L. Urban is used to estimate eta. It is a fit to
605// the results of L. Kim & al: Phys Rev. A33,3002 (1986)
606
607{
608 static const G4double K = 132.9416*eV;
609 static const G4double a1 = 4.15e-1, a3 = 2.10e-3, a5 = 54.0e-5;
610
611 G4double x = log(kineticEnergy/(K*Z*Z));
612 G4double x2 = x*x;
613 G4double x3 = x2*x;
614 G4double eta = 0.5 + atan(a1*x + a3*x3 + a5*x3*x2)/pi ;
615 G4double alfa = (1. - eta)/eta;
616 return eta*pow((1. - cut/kineticEnergy), alfa);
617}
618
619//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
620
621G4DataVector* G4eBremsstrahlungModel::ComputePartialSumSigma(
622 const G4Material* material,
623 G4double kineticEnergy,
624 G4double cut)
625
626// Build the table of cross section per element.
627//The table is built for MATERIALS.
628// This table is used by DoIt to select randomly an element in the material.
629{
630 G4int nElements = material->GetNumberOfElements();
631 const G4ElementVector* theElementVector = material->GetElementVector();
632 const G4double* theAtomNumDensityVector = material->GetAtomicNumDensityVector();
633 G4double dum = 0.;
634
635 G4DataVector* dv = new G4DataVector();
636
637 G4double cross = 0.0;
638
639 for (G4int i=0; i<nElements; i++ ) {
640
641 cross += theAtomNumDensityVector[i] * ComputeCrossSectionPerAtom( particle,
642 kineticEnergy, (*theElementVector)[i]->GetZ(), dum, cut);
643 dv->push_back(cross);
644 }
645 return dv;
646}
647
648//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
649
650void G4eBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
651 const G4MaterialCutsCouple* couple,
652 const G4DynamicParticle* dp,
653 G4double tmin,
654 G4double maxEnergy)
655// The emitted gamma energy is sampled using a parametrized formula
656// from L. Urban.
657// This parametrization is derived from :
658// cross-section values of Seltzer and Berger for electron energies
659// 1 keV - 10 GeV,
660// screened Bethe Heilter differential cross section above 10 GeV,
661// Migdal corrections in both case.
662// Seltzer & Berger: Nim B 12:95 (1985)
663// Nelson, Hirayama & Rogers: Technical report 265 SLAC (1985)
664// Migdal: Phys Rev 103:1811 (1956); Messel & Crawford: Pergamon Press (1970)
665//
666// A modified version of the random number techniques of Butcher&Messel is used
667// (Nuc Phys 20(1960),15).
668{
669 G4double kineticEnergy = dp->GetKineticEnergy();
670 G4double tmax = min(maxEnergy, kineticEnergy);
671 if(tmin >= tmax) return;
672
673//
674// GEANT4 internal units.
675//
676 static const G4double
677 ah10 = 4.67733E+00, ah11 =-6.19012E-01, ah12 = 2.02225E-02,
678 ah20 =-7.34101E+00, ah21 = 1.00462E+00, ah22 =-3.20985E-02,
679 ah30 = 2.93119E+00, ah31 =-4.03761E-01, ah32 = 1.25153E-02;
680
681 static const G4double
682 bh10 = 4.23071E+00, bh11 =-6.10995E-01, bh12 = 1.95531E-02,
683 bh20 =-7.12527E+00, bh21 = 9.69160E-01, bh22 =-2.74255E-02,
684 bh30 = 2.69925E+00, bh31 =-3.63283E-01, bh32 = 9.55316E-03;
685
686 static const G4double
687 al00 =-2.05398E+00, al01 = 2.38815E-02, al02 = 5.25483E-04,
688 al10 =-7.69748E-02, al11 =-6.91499E-02, al12 = 2.22453E-03,
689 al20 = 4.06463E-02, al21 =-1.01281E-02, al22 = 3.40919E-04;
690
691 static const G4double
692 bl00 = 1.04133E+00, bl01 =-9.43291E-03, bl02 =-4.54758E-04,
693 bl10 = 1.19253E-01, bl11 = 4.07467E-02, bl12 =-1.30718E-03,
694 bl20 =-1.59391E-02, bl21 = 7.27752E-03, bl22 =-1.94405E-04;
695
696 static const G4double tlow = 1.*MeV;
697
698 G4double gammaEnergy;
699 G4bool LPMOK = false;
700 const G4Material* material = couple->GetMaterial();
701
702 // select randomly one element constituing the material
703 const G4Element* anElement = SelectRandomAtom(couple);
704
705 // Extract Z factors for this Element
706 G4double lnZ = 3.*(anElement->GetIonisation()->GetlogZ3());
707 G4double FZ = lnZ* (4.- 0.55*lnZ);
708 G4double ZZ = anElement->GetIonisation()->GetZZ3();
709
710 // limits of the energy sampling
711 G4double totalEnergy = kineticEnergy + electron_mass_c2;
712 G4ThreeVector direction = dp->GetMomentumDirection();
713 G4double xmin = tmin/kineticEnergy;
714 G4double xmax = tmax/kineticEnergy;
715 G4double kappa = 0.0;
716 if(xmax >= 1.) xmax = 1.;
717 else kappa = log(xmax)/log(xmin);
718 G4double epsilmin = tmin/totalEnergy;
719 G4double epsilmax = tmax/totalEnergy;
720
721 // Migdal factor
722 G4double MigdalFactor = (material->GetElectronDensity())*MigdalConstant
723 / (epsilmax*epsilmax);
724
725 G4double x, epsil, greject, migdal, grejmax, q;
726 G4double U = log(kineticEnergy/electron_mass_c2);
727 G4double U2 = U*U;
728
729 // precalculated parameters
730 G4double ah, bh;
731 G4double screenfac = 0.0;
732
733 if (kineticEnergy > tlow) {
734
735 G4double ah1 = ah10 + ZZ* (ah11 + ZZ* ah12);
736 G4double ah2 = ah20 + ZZ* (ah21 + ZZ* ah22);
737 G4double ah3 = ah30 + ZZ* (ah31 + ZZ* ah32);
738
739 G4double bh1 = bh10 + ZZ* (bh11 + ZZ* bh12);
740 G4double bh2 = bh20 + ZZ* (bh21 + ZZ* bh22);
741 G4double bh3 = bh30 + ZZ* (bh31 + ZZ* bh32);
742
743 ah = 1. + (ah1*U2 + ah2*U + ah3) / (U2*U);
744 bh = 0.75 + (bh1*U2 + bh2*U + bh3) / (U2*U);
745
746 // limit of the screening variable
747 screenfac =
748 136.*electron_mass_c2/((anElement->GetIonisation()->GetZ3())*totalEnergy);
749 G4double screenmin = screenfac*epsilmin/(1.-epsilmin);
750
751 // Compute the maximum of the rejection function
752 G4double F1 = max(ScreenFunction1(screenmin) - FZ ,0.);
753 G4double F2 = max(ScreenFunction2(screenmin) - FZ ,0.);
754 grejmax = (F1 - epsilmin* (F1*ah - bh*epsilmin*F2))/(42.392 - FZ);
755
756 } else {
757
758 G4double al0 = al00 + ZZ* (al01 + ZZ* al02);
759 G4double al1 = al10 + ZZ* (al11 + ZZ* al12);
760 G4double al2 = al20 + ZZ* (al21 + ZZ* al22);
761
762 G4double bl0 = bl00 + ZZ* (bl01 + ZZ* bl02);
763 G4double bl1 = bl10 + ZZ* (bl11 + ZZ* bl12);
764 G4double bl2 = bl20 + ZZ* (bl21 + ZZ* bl22);
765
766 ah = al0 + al1*U + al2*U2;
767 bh = bl0 + bl1*U + bl2*U2;
768
769 // Compute the maximum of the rejection function
770 grejmax = max(1. + xmin* (ah + bh*xmin), 1.+ah+bh);
771 G4double xm = -ah/(2.*bh);
772 if ( xmin < xm && xm < xmax) grejmax = max(grejmax, 1.+ xm* (ah + bh*xm));
773 }
774
775 //
776 // sample the energy rate of the emitted gamma for e- kin energy > 1 MeV
777 //
778
779 do {
780 if (kineticEnergy > tlow) {
781 do {
782 q = G4UniformRand();
783 x = pow(xmin, q + kappa*(1.0 - q));
784 epsil = x*kineticEnergy/totalEnergy;
785 G4double screenvar = screenfac*epsil/(1.0-epsil);
786 G4double F1 = max(ScreenFunction1(screenvar) - FZ ,0.);
787 G4double F2 = max(ScreenFunction2(screenvar) - FZ ,0.);
788 migdal = (1. + MigdalFactor)/(1. + MigdalFactor/(x*x));
789 greject = migdal*(F1 - epsil* (ah*F1 - bh*epsil*F2))/(42.392 - FZ);
790 /*
791 if ( greject > grejmax ) {
792 G4cout << "### G4eBremsstrahlungModel Warning: Majoranta exceeded! "
793 << greject << " > " << grejmax
794 << " x= " << x
795 << " e= " << kineticEnergy
796 << G4endl;
797 }
798 */
799 } while( greject < G4UniformRand()*grejmax );
800
801 } else {
802
803 do {
804 q = G4UniformRand();
805 x = pow(xmin, q + kappa*(1.0 - q));
806 migdal = (1. + MigdalFactor)/(1. + MigdalFactor/(x*x));
807 greject = migdal*(1. + x* (ah + bh*x));
808 /*
809 if ( greject > grejmax ) {
810 G4cout << "### G4eBremsstrahlungModel Warning: Majoranta exceeded! "
811 << greject << " > " << grejmax
812 << " x= " << x
813 << " e= " << kineticEnergy
814 << G4endl;
815 }
816 */
817 } while( greject < G4UniformRand()*grejmax );
818 }
819 /*
820 if(x > 0.999) {
821 G4cout << "### G4eBremsstrahlungModel Warning: e= " << kineticEnergy
822 << " tlow= " << tlow
823 << " x= " << x
824 << " greject= " << greject
825 << " grejmax= " << grejmax
826 << " migdal= " << migdal
827 << G4endl;
828 // if(x >= 1.0) G4Exception("X=1");
829 }
830 */
831 gammaEnergy = x*kineticEnergy;
832
833 if (LPMFlag()) {
834 // take into account the supression due to the LPM effect
835 if (G4UniformRand() <= SupressionFunction(material,kineticEnergy,
836 gammaEnergy))
837 LPMOK = true;
838 }
839 else LPMOK = true;
840
841 } while (!LPMOK);
842
843 //
844 // angles of the emitted gamma. ( Z - axis along the parent particle)
845 //
846 // universal distribution suggested by L. Urban
847 // (Geant3 manual (1993) Phys211),
848 // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
849
850 G4double u;
851 const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
852
853 if (9./(9.+d) > G4UniformRand()) u = - log(G4UniformRand()*G4UniformRand())/a1;
854 else u = - log(G4UniformRand()*G4UniformRand())/a2;
855
856 G4double theta = u*electron_mass_c2/totalEnergy;
857
858 G4double sint = sin(theta);
859
860 G4double phi = twopi * G4UniformRand() ;
861
862 G4ThreeVector gammaDirection(sint*cos(phi),sint*sin(phi), cos(theta));
863 gammaDirection.rotateUz(direction);
864
865 // create G4DynamicParticle object for the Gamma
866 G4DynamicParticle* g = new G4DynamicParticle(theGamma,gammaDirection,
867 gammaEnergy);
868 vdp->push_back(g);
869
870 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
871 G4ThreeVector dir = totMomentum*direction - gammaEnergy*gammaDirection;
872 direction = dir.unit();
873
874 // energy of primary
875 G4double finalE = kineticEnergy - gammaEnergy;
876
877 // stop tracking and create new secondary instead of primary
878 if(gammaEnergy > SecondaryThreshold()) {
879 fParticleChange->ProposeTrackStatus(fStopAndKill);
880 fParticleChange->SetProposedKineticEnergy(0.0);
881 G4DynamicParticle* el =
882 new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
883 direction, finalE);
884 vdp->push_back(el);
885
886 // continue tracking
887 } else {
888 fParticleChange->SetProposedMomentumDirection(direction);
889 fParticleChange->SetProposedKineticEnergy(finalE);
890 }
891}
892
893//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
894
895const G4Element* G4eBremsstrahlungModel::SelectRandomAtom(
896 const G4MaterialCutsCouple* couple)
897{
898 // select randomly 1 element within the material
899
900 const G4Material* material = couple->GetMaterial();
901 G4int nElements = material->GetNumberOfElements();
902 const G4ElementVector* theElementVector = material->GetElementVector();
903
904 const G4Element* elm = 0;
905
906 if(1 < nElements) {
907
908 G4DataVector* dv = partialSumSigma[couple->GetIndex()];
909 G4double rval = G4UniformRand()*((*dv)[nElements-1]);
910
911 for (G4int i=0; i<nElements; i++) {
912 if (rval <= (*dv)[i]) elm = (*theElementVector)[i];
913 }
914 if(!elm) {
915 G4cout << "G4eBremsstrahlungModel::SelectRandomAtom: Warning -"
916 << " no elements found in "
917 << material->GetName()
918 << G4endl;
919 elm = (*theElementVector)[0];
920 }
921 } else elm = (*theElementVector)[0];
922
923 SetCurrentElement(elm);
924 return elm;
925}
926
927//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
928
929G4double G4eBremsstrahlungModel::SupressionFunction(const G4Material* material,
930 G4double kineticEnergy, G4double gammaEnergy)
931{
932 // supression due to the LPM effect+polarisation of the medium/
933 // supression due to the polarisation alone
934
935
936 G4double totEnergy = kineticEnergy+electron_mass_c2 ;
937 G4double totEnergySquare = totEnergy*totEnergy ;
938
939 G4double LPMEnergy = LPMconstant*(material->GetRadlen()) ;
940
941 G4double gammaEnergySquare = gammaEnergy*gammaEnergy ;
942
943 G4double electronDensity = material->GetElectronDensity();
944
945 G4double sp = gammaEnergySquare/
946 (gammaEnergySquare+MigdalConstant*totEnergySquare*electronDensity);
947
948 G4double supr = 1.0;
949
950 if (LPMFlag()) {
951
952 G4double s2lpm = LPMEnergy*gammaEnergy/totEnergySquare;
953
954 if (s2lpm < 1.) {
955
956 G4double LPMgEnergyLimit = totEnergySquare/LPMEnergy ;
957 G4double LPMgEnergyLimit2 = LPMgEnergyLimit*LPMgEnergyLimit;
958 G4double splim = LPMgEnergyLimit2/
959 (LPMgEnergyLimit2+MigdalConstant*totEnergySquare*electronDensity);
960 G4double w = 1.+1./splim ;
961
962 if ((1.-sp) < 1.e-6) w = s2lpm*(3.-sp);
963 else w = s2lpm*(1.+1./sp);
964
965 supr = (sqrt(w*w+4.*s2lpm)-w)/(sqrt(w*w+4.)-w) ;
966 supr /= sp;
967 }
968
969 }
970 return supr;
971}
972
973//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
974
975
Note: See TracBrowser for help on using the repository browser.