source: trunk/source/processes/electromagnetic/standard/src/G4eBremsstrahlungRelModel.cc@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 17.5 KB
RevLine 
[968]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//
[1055]26// $Id: G4eBremsstrahlungRelModel.cc,v 1.14 2009/04/09 18:41:18 vnivanch Exp $
[1196]27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[968]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4eBremsstrahlungRelModel
35//
36// Author: Andreas Schaelicke
37//
38// Creation date: 12.08.2008
39//
40// Modifications:
41//
42// 13.11.08 add SetLPMflag and SetLPMconstant methods
43// 13.11.08 change default LPMconstant value
44//
45// Main References:
46// Y.-S.Tsai, Rev. Mod. Phys. 46 (1974) 815; Rev. Mod. Phys. 49 (1977) 421.
47// S.Klein, Rev. Mod. Phys. 71 (1999) 1501.
48// T.Stanev et.al., Phys. Rev. D25 (1982) 1291.
49// M.L.Ter-Mikaelian, High-energy Electromagnetic Processes in Condensed Media, Wiley, 1972.
50//
51// -------------------------------------------------------------------
52//
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
56#include "G4eBremsstrahlungRelModel.hh"
57#include "G4Electron.hh"
58#include "G4Positron.hh"
59#include "G4Gamma.hh"
60#include "Randomize.hh"
61#include "G4Material.hh"
62#include "G4Element.hh"
63#include "G4ElementVector.hh"
64#include "G4ProductionCutsTable.hh"
65#include "G4ParticleChangeForLoss.hh"
66#include "G4LossTableManager.hh"
67
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
71const G4double G4eBremsstrahlungRelModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
72 0.5917, 0.7628, 0.8983, 0.9801 };
73const G4double G4eBremsstrahlungRelModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
74 0.1813, 0.1569, 0.1112, 0.0506 };
75const G4double G4eBremsstrahlungRelModel::Fel_light[] = {0., 5.31 , 4.79 , 4.74 , 4.71} ;
76const G4double G4eBremsstrahlungRelModel::Finel_light[] = {0., 6.144 , 5.621 , 5.805 , 5.924} ;
77
78
79using namespace std;
80
81G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel(const G4ParticleDefinition* p,
82 const G4String& name)
83 : G4VEmModel(name),
84 particle(0),
85 fXiLPM(0), fPhiLPM(0), fGLPM(0),
86 isElectron(true),
87 fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
88 fLPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)*0.5),
89 bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
90 use_completescreening(true),isInitialised(false)
91{
92 if(p) SetParticle(p);
93 theGamma = G4Gamma::Gamma();
94
95 minThreshold = 1.0*keV;
96 SetLowEnergyLimit(GeV);
97
98 nist = G4NistManager::Instance();
99 InitialiseConstants();
100
101 SetLPMFlag(true);
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
106void G4eBremsstrahlungRelModel::InitialiseConstants()
107{
108 facFel = log(184.15);
109 facFinel = log(1194.);
110
111 preS1 = 1./(184.15*184.15);
112 logTwo = log(2.);
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
117G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel()
118{
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
123void G4eBremsstrahlungRelModel::SetParticle(const G4ParticleDefinition* p)
124{
125 particle = p;
126 particleMass = p->GetPDGMass();
127 if(p == G4Electron::Electron()) isElectron = true;
128 else isElectron = false;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
133G4double G4eBremsstrahlungRelModel::MinEnergyCut(const G4ParticleDefinition*,
134 const G4MaterialCutsCouple*)
135{
136 return minThreshold;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
141void G4eBremsstrahlungRelModel::SetupForMaterial(const G4ParticleDefinition*,
142 const G4Material* mat, G4double kineticEnergy)
143{
144 densityFactor = mat->GetElectronDensity()*fMigdalConstant;
145 lpmEnergy = mat->GetRadlen()*fLPMconstant;
146
147 // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
148 if (LPMFlag())
149 energyThresholdLPM=sqrt(densityFactor)*lpmEnergy;
150 else
151 energyThresholdLPM=1.e39; // i.e. do not use LPM effect
152
153 // calculate threshold for density effect
154 kinEnergy = kineticEnergy;
155 totalEnergy = kineticEnergy + particleMass;
156 densityCorr = densityFactor*totalEnergy*totalEnergy;
157
158 // define critical gamma energies (important for integration/dicing)
159 klpm=totalEnergy*totalEnergy/lpmEnergy;
160 kp=sqrt(densityCorr);
161
162}
163
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167void G4eBremsstrahlungRelModel::Initialise(const G4ParticleDefinition* p,
168 const G4DataVector& cuts)
169{
170 if(p) SetParticle(p);
171
172 highKinEnergy = HighEnergyLimit();
173 lowKinEnergy = LowEnergyLimit();
174
175 currentZ = 0.;
176
177 InitialiseElementSelectors(p, cuts);
178
179 if(isInitialised) return;
[1055]180 fParticleChange = GetParticleChangeForLoss();
[968]181 isInitialised = true;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
186G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(
187 const G4Material* material,
188 const G4ParticleDefinition* p,
189 G4double kineticEnergy,
190 G4double cutEnergy)
191{
192 if(!particle) SetParticle(p);
193 if(kineticEnergy < lowKinEnergy) return 0.0;
194 G4double cut = std::min(cutEnergy, kineticEnergy);
195 if(cut == 0.0) return 0.0;
196
197 SetupForMaterial(particle, material,kineticEnergy);
198
199 const G4ElementVector* theElementVector = material->GetElementVector();
200 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
201
202 G4double dedx = 0.0;
203
204 // loop for elements in the material
205 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
206
207 G4VEmModel::SetCurrentElement((*theElementVector)[i]);
208 SetCurrentElement((*theElementVector)[i]->GetZ());
209
210 dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
211 }
212 dedx *= bremFactor;
213
214
215 return dedx;
216}
217
218//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
219
220G4double G4eBremsstrahlungRelModel::ComputeBremLoss(G4double cut)
221{
222 G4double loss = 0.0;
223
224 // number of intervals and integration step
225 G4double vcut = cut/totalEnergy;
226 G4int n = (G4int)(20*vcut) + 3;
227 G4double delta = vcut/G4double(n);
228
229 G4double e0 = 0.0;
230 G4double xs;
231
232 // integration
233 for(G4int l=0; l<n; l++) {
234
235 for(G4int i=0; i<8; i++) {
236
237 G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
238
239 if(totalEnergy > energyThresholdLPM) {
240 xs = ComputeRelDXSectionPerAtom(eg);
241 } else {
242 xs = ComputeDXSectionPerAtom(eg);
243 }
244 loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
245 }
246 e0 += delta;
247 }
248
249 loss *= delta*totalEnergy;
250
251 return loss;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
255
256G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom(
257 const G4ParticleDefinition* p,
258 G4double kineticEnergy,
259 G4double Z, G4double,
260 G4double cutEnergy,
261 G4double maxEnergy)
262{
263 if(!particle) SetParticle(p);
264 if(kineticEnergy < lowKinEnergy) return 0.0;
265 G4double cut = std::min(cutEnergy, kineticEnergy);
266 G4double tmax = std::min(maxEnergy, kineticEnergy);
267
268 if(cut >= tmax) return 0.0;
269
270 SetCurrentElement(Z);
271
272 G4double cross = ComputeXSectionPerAtom(cut);
273
274 // allow partial integration
275 if(tmax < kinEnergy) cross -= ComputeXSectionPerAtom(tmax);
276
277 cross *= Z*Z*bremFactor;
278
279 return cross;
280}
281
282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
283
284
285G4double G4eBremsstrahlungRelModel::ComputeXSectionPerAtom(G4double cut)
286{
287 G4double cross = 0.0;
288
289 // number of intervals and integration step
290 G4double vcut = log(cut/totalEnergy);
291 G4double vmax = log(kinEnergy/totalEnergy);
292 G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
293 // n=1; // integration test
294 G4double delta = (vmax - vcut)/G4double(n);
295
296 G4double e0 = vcut;
297 G4double xs;
298
299 // integration
300 for(G4int l=0; l<n; l++) {
301
302 for(G4int i=0; i<8; i++) {
303
304 G4double eg = exp(e0 + xgi[i]*delta)*totalEnergy;
305
306 if(totalEnergy > energyThresholdLPM) {
307 xs = ComputeRelDXSectionPerAtom(eg);
308 } else {
309 xs = ComputeDXSectionPerAtom(eg);
310 }
311 cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
312 }
313 e0 += delta;
314 }
315
316 cross *= delta;
317
318 return cross;
319}
320
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
322void G4eBremsstrahlungRelModel::CalcLPMFunctions(G4double k)
323{
324 // *** calculate lpm variable s & sprime ***
325 // Klein eqs. (78) & (79)
326 G4double sprime = sqrt(0.125*k*lpmEnergy/(totalEnergy*(totalEnergy-k)));
327
328 G4double s1 = preS1*z23;
329 G4double logS1 = 2./3.*lnZ-2.*facFel;
330 G4double logTS1 = logTwo+logS1;
331
332 xiLPM = 2.;
333
334 if (sprime>1)
335 xiLPM = 1.;
336 else if (sprime>sqrt(2.)*s1) {
337 G4double h = log(sprime)/logTS1;
338 xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
339 }
340
341 G4double s = sprime/sqrt(xiLPM);
342
343 // *** merging with density effect*** should be only necessary in region "close to" kp, e.g. k<100*kp
344 // using Ter-Mikaelian eq. (20.9)
345 G4double k2 = k*k;
346 s = s * (1 + (densityCorr/k2) );
347
348 // recalculate Xi using modified s above
349 // Klein eq. (75)
350 xiLPM = 1.;
351 if (s<=s1) xiLPM = 2.;
352 else if ( (s1<s) && (s<=1) ) xiLPM = 1. + log(s)/logS1;
353
354
355 // *** calculate supression functions phi and G ***
356 // Klein eqs. (77)
357 G4double s2=s*s;
358 G4double s3=s*s2;
359 G4double s4=s2*s2;
360
361 if (s<0.1) {
362 // high suppression limit
363 phiLPM = 6.*s - 18.84955592153876*s2 + 39.47841760435743*s3
364 - 57.69873135166053*s4;
365 gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
366 }
367 else if (s<1.9516) {
368 // intermediate suppression
369 // using eq.77 approxim. valid s<2.
370 phiLPM = 1.-exp(-6.*s*(1.+(3.-pi)*s)
371 +s3/(0.623+0.795*s+0.658*s2));
372 if (s<0.415827397755) {
373 // using eq.77 approxim. valid 0.07<s<2
374 G4double psiLPM = 1-exp(-4*s-8*s2/(1+3.936*s+4.97*s2-0.05*s3+7.50*s4));
375 gLPM = 3*psiLPM-2*phiLPM;
376 }
377 else {
378 // using alternative parametrisiation
379 G4double pre = -0.16072300849123999 + s*3.7550300067531581 + s2*-1.7981383069010097
380 + s3*0.67282686077812381 + s4*-0.1207722909879257;
381 gLPM = tanh(pre);
382 }
383 }
384 else {
385 // low suppression limit valid s>2.
386 phiLPM = 1. - 0.0119048/s4;
387 gLPM = 1. - 0.0230655/s4;
388 }
389
390 // *** make sure suppression is smaller than 1 ***
391 // *** caused by Migdal approximation in xi ***
392 if (xiLPM*phiLPM>1. || s>0.57) xiLPM=1./phiLPM;
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
396
397
398G4double G4eBremsstrahlungRelModel::ComputeRelDXSectionPerAtom(G4double gammaEnergy)
399// Ultra relativistic model
400// only valid for very high energies, but includes LPM suppression
401// * complete screening
402{
403 if(gammaEnergy < 0.0) return 0.0;
404
405 G4double y = gammaEnergy/totalEnergy;
406 G4double y2 = y*y*.25;
407 G4double yone2 = (1.-y+2.*y2);
408
409 // ** form factors complete screening case **
410
411 // ** calc LPM functions -- include ter-mikaelian merging with density effect **
412 // G4double xiLPM, gLPM, phiLPM; // to be made member variables !!!
413 CalcLPMFunctions(gammaEnergy);
414
415 G4double mainLPM = xiLPM*(y2 * gLPM + yone2*phiLPM) * ( (Fel-fCoulomb) + Finel/currentZ );
416 G4double secondTerm = (1.-y)/12.*(1.+1./currentZ);
417
418 G4double cross = mainLPM+secondTerm;
419 return cross;
420}
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
423
424G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom(G4double gammaEnergy)
425// Relativistic model
426// only valid for high energies (and if LPM suppression does not play a role)
427// * screening according to thomas-fermi-Model (only valid for Z>5)
428// * no LPM effect
429{
430
431 if(gammaEnergy < 0.0) return 0.0;
432
433 G4double y = gammaEnergy/totalEnergy;
434
435 G4double main=0.,secondTerm=0.;
436
437 if (use_completescreening|| currentZ<5) {
438 // ** form factors complete screening case **
439 main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
440 secondTerm = (1.-y)/12.*(1.+1./currentZ);
441 }
442 else {
443 // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
444 G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
445 G4double gg=dd*z13;
446 G4double eps=dd*z23;
447 G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ);
448 G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ);
449
450 main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
451 secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
452 }
453 G4double cross = main+secondTerm;
454 return cross;
455}
456
457//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458
459void G4eBremsstrahlungRelModel::SampleSecondaries(
460 std::vector<G4DynamicParticle*>* vdp,
461 const G4MaterialCutsCouple* couple,
462 const G4DynamicParticle* dp,
463 G4double cutEnergy,
464 G4double maxEnergy)
465{
466 G4double kineticEnergy = dp->GetKineticEnergy();
467 if(kineticEnergy < lowKinEnergy) return;
468 G4double cut = std::min(cutEnergy, kineticEnergy);
469 G4double emax = std::min(maxEnergy, kineticEnergy);
470 if(cut >= emax) return;
471
472 SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
473
474 const G4Element* elm =
475 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
476 SetCurrentElement(elm->GetZ());
477
478 kinEnergy = kineticEnergy;
479 totalEnergy = kineticEnergy + particleMass;
480 densityCorr = densityFactor*totalEnergy*totalEnergy;
481 G4ThreeVector direction = dp->GetMomentumDirection();
482
483 // G4double fmax= fMax;
484 G4bool highe = true;
485 if(totalEnergy < energyThresholdLPM) highe = false;
486
487 G4double xmin = log(cut*cut + densityCorr);
488 G4double xmax = log(emax*emax + densityCorr);
489 G4double gammaEnergy, f, x;
490
491 do {
492 x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
493 if(x < 0.0) x = 0.0;
494 gammaEnergy = sqrt(x);
495 if(highe) f = ComputeRelDXSectionPerAtom(gammaEnergy);
496 else f = ComputeDXSectionPerAtom(gammaEnergy);
497
498 if ( f > fMax ) {
499 G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
500 << f << " > " << fMax
501 << " Egamma(MeV)= " << gammaEnergy
502 << " E(mEV)= " << kineticEnergy
503 << G4endl;
504 }
505
506 } while (f < fMax*G4UniformRand());
507
508 //
509 // angles of the emitted gamma. ( Z - axis along the parent particle)
510 //
511 // universal distribution suggested by L. Urban
512 // (Geant3 manual (1993) Phys211),
513 // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
514
515 G4double u;
516 const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
517
518 if (9./(9.+d) > G4UniformRand()) u = - log(G4UniformRand()*G4UniformRand())/a1;
519 else u = - log(G4UniformRand()*G4UniformRand())/a2;
520
521 G4double theta = u*particleMass/totalEnergy;
522 G4double sint = sin(theta);
523 G4double phi = twopi * G4UniformRand();
524 G4ThreeVector gammaDirection(sint*cos(phi),sint*sin(phi), cos(theta));
525 gammaDirection.rotateUz(direction);
526
527 // create G4DynamicParticle object for the Gamma
528 G4DynamicParticle* g = new G4DynamicParticle(theGamma,gammaDirection,
529 gammaEnergy);
530 vdp->push_back(g);
531
532 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
533 G4ThreeVector dir = totMomentum*direction - gammaEnergy*gammaDirection;
534 direction = dir.unit();
535
536 // energy of primary
537 G4double finalE = kineticEnergy - gammaEnergy;
538
539 // stop tracking and create new secondary instead of primary
540 if(gammaEnergy > SecondaryThreshold()) {
541 fParticleChange->ProposeTrackStatus(fStopAndKill);
542 fParticleChange->SetProposedKineticEnergy(0.0);
543 G4DynamicParticle* el =
544 new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
545 direction, finalE);
546 vdp->push_back(el);
547
548 // continue tracking
549 } else {
550 fParticleChange->SetProposedMomentumDirection(direction);
551 fParticleChange->SetProposedKineticEnergy(finalE);
552 }
553}
554
555//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
556
557
Note: See TracBrowser for help on using the repository browser.