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

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

fichier ajoutes

File size: 17.6 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: G4eBremsstrahlungRelModel.cc,v 1.13 2009/02/20 12:06:37 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
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;
180
181  if(pParticleChange) {
182    fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
183  } else {
184    fParticleChange = new G4ParticleChangeForLoss();
185  }
186  isInitialised = true;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190
191G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(
192                                             const G4Material* material,
193                                             const G4ParticleDefinition* p,
194                                                   G4double kineticEnergy,
195                                                   G4double cutEnergy)
196{
197  if(!particle) SetParticle(p);
198  if(kineticEnergy < lowKinEnergy) return 0.0;
199  G4double cut = std::min(cutEnergy, kineticEnergy);
200  if(cut == 0.0) return 0.0;
201
202  SetupForMaterial(particle, material,kineticEnergy);
203
204  const G4ElementVector* theElementVector = material->GetElementVector();
205  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
206
207  G4double dedx = 0.0;
208
209  //  loop for elements in the material
210  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
211
212    G4VEmModel::SetCurrentElement((*theElementVector)[i]);
213    SetCurrentElement((*theElementVector)[i]->GetZ());
214
215    dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
216  }
217  dedx *= bremFactor;
218
219
220  return dedx;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
225G4double G4eBremsstrahlungRelModel::ComputeBremLoss(G4double cut)
226{
227  G4double loss = 0.0;
228
229  // number of intervals and integration step
230  G4double vcut = cut/totalEnergy;
231  G4int n = (G4int)(20*vcut) + 3;
232  G4double delta = vcut/G4double(n);
233
234  G4double e0 = 0.0;
235  G4double xs; 
236
237  // integration
238  for(G4int l=0; l<n; l++) {
239
240    for(G4int i=0; i<8; i++) {
241
242      G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
243
244      if(totalEnergy > energyThresholdLPM) {
245        xs = ComputeRelDXSectionPerAtom(eg);
246      } else {
247        xs = ComputeDXSectionPerAtom(eg);
248      }
249      loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
250    }
251    e0 += delta;
252  }
253
254  loss *= delta*totalEnergy;
255
256  return loss;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260
261G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom(
262                                              const G4ParticleDefinition* p,
263                                              G4double kineticEnergy, 
264                                              G4double Z,   G4double,
265                                              G4double cutEnergy, 
266                                              G4double maxEnergy)
267{
268  if(!particle) SetParticle(p);
269  if(kineticEnergy < lowKinEnergy) return 0.0;
270  G4double cut  = std::min(cutEnergy, kineticEnergy);
271  G4double tmax = std::min(maxEnergy, kineticEnergy);
272
273  if(cut >= tmax) return 0.0;
274
275  SetCurrentElement(Z);
276
277  G4double cross = ComputeXSectionPerAtom(cut);
278
279  // allow partial integration
280  if(tmax < kinEnergy) cross -= ComputeXSectionPerAtom(tmax);
281 
282  cross *= Z*Z*bremFactor;
283
284  return cross;
285}
286 
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
288
289
290G4double G4eBremsstrahlungRelModel::ComputeXSectionPerAtom(G4double cut)
291{
292  G4double cross = 0.0;
293
294  // number of intervals and integration step
295  G4double vcut = log(cut/totalEnergy);
296  G4double vmax = log(kinEnergy/totalEnergy);
297  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
298  //  n=1; //  integration test
299  G4double delta = (vmax - vcut)/G4double(n);
300
301  G4double e0 = vcut;
302  G4double xs; 
303
304  // integration
305  for(G4int l=0; l<n; l++) {
306
307    for(G4int i=0; i<8; i++) {
308
309      G4double eg = exp(e0 + xgi[i]*delta)*totalEnergy;
310
311      if(totalEnergy > energyThresholdLPM) {
312        xs = ComputeRelDXSectionPerAtom(eg);
313      } else {
314        xs = ComputeDXSectionPerAtom(eg);
315      }
316      cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
317    }
318    e0 += delta;
319  }
320
321  cross *= delta;
322
323  return cross;
324}
325
326//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
327void  G4eBremsstrahlungRelModel::CalcLPMFunctions(G4double k)
328{
329  // *** calculate lpm variable s & sprime ***
330  // Klein eqs. (78) & (79)
331  G4double sprime = sqrt(0.125*k*lpmEnergy/(totalEnergy*(totalEnergy-k)));
332
333  G4double s1 = preS1*z23;
334  G4double logS1 = 2./3.*lnZ-2.*facFel;
335  G4double logTS1 = logTwo+logS1;
336
337  xiLPM = 2.;
338
339  if (sprime>1) 
340    xiLPM = 1.;
341  else if (sprime>sqrt(2.)*s1) {
342    G4double h  = log(sprime)/logTS1;
343    xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
344  }
345
346  G4double s = sprime/sqrt(xiLPM); 
347
348  // *** merging with density effect***  should be only necessary in region "close to" kp, e.g. k<100*kp
349  // using Ter-Mikaelian eq. (20.9)
350  G4double k2 = k*k;
351  s = s * (1 + (densityCorr/k2) );
352
353  // recalculate Xi using modified s above
354  // Klein eq. (75)
355  xiLPM = 1.;
356  if (s<=s1) xiLPM = 2.;
357  else if ( (s1<s) && (s<=1) ) xiLPM = 1. + log(s)/logS1;
358 
359
360  // *** calculate supression functions phi and G ***
361  // Klein eqs. (77)
362  G4double s2=s*s;
363  G4double s3=s*s2;
364  G4double s4=s2*s2;
365
366  if (s<0.1) {
367    // high suppression limit
368    phiLPM = 6.*s - 18.84955592153876*s2 + 39.47841760435743*s3
369      - 57.69873135166053*s4;
370    gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
371  }
372  else if (s<1.9516) {
373    // intermediate suppression
374    // using eq.77 approxim. valid s<2.     
375    phiLPM = 1.-exp(-6.*s*(1.+(3.-pi)*s)
376                +s3/(0.623+0.795*s+0.658*s2));
377    if (s<0.415827397755) {
378      // using eq.77 approxim. valid 0.07<s<2
379      G4double psiLPM = 1-exp(-4*s-8*s2/(1+3.936*s+4.97*s2-0.05*s3+7.50*s4));
380      gLPM = 3*psiLPM-2*phiLPM;
381    }
382    else {
383      // using alternative parametrisiation
384      G4double pre = -0.16072300849123999 + s*3.7550300067531581 + s2*-1.7981383069010097 
385        + s3*0.67282686077812381 + s4*-0.1207722909879257;
386      gLPM = tanh(pre);
387    }
388  }
389  else {
390    // low suppression limit valid s>2.
391    phiLPM = 1. - 0.0119048/s4;
392    gLPM = 1. - 0.0230655/s4;
393  }
394
395  // *** make sure suppression is smaller than 1 ***
396  // *** caused by Migdal approximation in xi    ***
397  if (xiLPM*phiLPM>1. || s>0.57)  xiLPM=1./phiLPM;
398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401
402
403G4double G4eBremsstrahlungRelModel::ComputeRelDXSectionPerAtom(G4double gammaEnergy)
404// Ultra relativistic model
405//   only valid for very high energies, but includes LPM suppression
406//    * complete screening
407{
408  if(gammaEnergy < 0.0) return 0.0;
409
410  G4double y = gammaEnergy/totalEnergy;
411  G4double y2 = y*y*.25;
412  G4double yone2 = (1.-y+2.*y2);
413
414  // ** form factors complete screening case **     
415
416  // ** calc LPM functions -- include ter-mikaelian merging with density effect **
417  //  G4double xiLPM, gLPM, phiLPM;  // to be made member variables !!!
418  CalcLPMFunctions(gammaEnergy);
419
420  G4double mainLPM   = xiLPM*(y2 * gLPM + yone2*phiLPM) * ( (Fel-fCoulomb) + Finel/currentZ );
421  G4double secondTerm = (1.-y)/12.*(1.+1./currentZ);
422
423  G4double cross = mainLPM+secondTerm;
424  return cross;
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
429G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom(G4double gammaEnergy)
430// Relativistic model
431//  only valid for high energies (and if LPM suppression does not play a role)
432//  * screening according to thomas-fermi-Model (only valid for Z>5)
433//  * no LPM effect
434{
435
436  if(gammaEnergy < 0.0) return 0.0;
437
438  G4double y = gammaEnergy/totalEnergy;
439
440  G4double main=0.,secondTerm=0.;
441
442  if (use_completescreening|| currentZ<5) {
443    // ** form factors complete screening case **     
444    main   = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
445    secondTerm = (1.-y)/12.*(1.+1./currentZ);
446  }
447  else {
448    // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
449    G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
450    G4double gg=dd*z13;
451    G4double eps=dd*z23;
452    G4double phi1=Phi1(gg,currentZ),  phi1m2=Phi1M2(gg,currentZ);
453    G4double psi1=Psi1(eps,currentZ),  psi1m2=Psi1M2(eps,currentZ);
454   
455    main   = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
456    secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
457  }
458  G4double cross = main+secondTerm;
459  return cross;
460}
461
462//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
463
464void G4eBremsstrahlungRelModel::SampleSecondaries(
465                                      std::vector<G4DynamicParticle*>* vdp, 
466                                      const G4MaterialCutsCouple* couple,
467                                      const G4DynamicParticle* dp,
468                                      G4double cutEnergy,
469                                      G4double maxEnergy)
470{
471  G4double kineticEnergy = dp->GetKineticEnergy();
472  if(kineticEnergy < lowKinEnergy) return;
473  G4double cut  = std::min(cutEnergy, kineticEnergy);
474  G4double emax = std::min(maxEnergy, kineticEnergy);
475  if(cut >= emax) return;
476
477  SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
478
479  const G4Element* elm = 
480    SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
481  SetCurrentElement(elm->GetZ());
482
483  kinEnergy   = kineticEnergy;
484  totalEnergy = kineticEnergy + particleMass;
485  densityCorr = densityFactor*totalEnergy*totalEnergy;
486  G4ThreeVector direction = dp->GetMomentumDirection();
487
488  //  G4double fmax= fMax;
489  G4bool highe = true;
490  if(totalEnergy < energyThresholdLPM) highe = false;
491 
492  G4double xmin = log(cut*cut + densityCorr);
493  G4double xmax = log(emax*emax  + densityCorr);
494  G4double gammaEnergy, f, x; 
495
496  do {
497    x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
498    if(x < 0.0) x = 0.0;
499    gammaEnergy = sqrt(x);
500    if(highe) f = ComputeRelDXSectionPerAtom(gammaEnergy);
501    else      f = ComputeDXSectionPerAtom(gammaEnergy);
502
503    if ( f > fMax ) {
504      G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
505             << f << " > " << fMax
506             << " Egamma(MeV)= " << gammaEnergy
507             << " E(mEV)= " << kineticEnergy
508             << G4endl;
509    }
510
511  } while (f < fMax*G4UniformRand());
512
513  //
514  //  angles of the emitted gamma. ( Z - axis along the parent particle)
515  //
516  //  universal distribution suggested by L. Urban
517  // (Geant3 manual (1993) Phys211),
518  //  derived from Tsai distribution (Rev Mod Phys 49,421(1977))
519
520  G4double u;
521  const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
522
523  if (9./(9.+d) > G4UniformRand()) u = - log(G4UniformRand()*G4UniformRand())/a1;
524     else                          u = - log(G4UniformRand()*G4UniformRand())/a2;
525
526  G4double theta = u*particleMass/totalEnergy;
527  G4double sint = sin(theta);
528  G4double phi = twopi * G4UniformRand();
529  G4ThreeVector gammaDirection(sint*cos(phi),sint*sin(phi), cos(theta));
530  gammaDirection.rotateUz(direction);
531
532  // create G4DynamicParticle object for the Gamma
533  G4DynamicParticle* g = new G4DynamicParticle(theGamma,gammaDirection,
534                                                        gammaEnergy);
535  vdp->push_back(g);
536 
537  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
538  G4ThreeVector dir = totMomentum*direction - gammaEnergy*gammaDirection;
539  direction = dir.unit();
540
541  // energy of primary
542  G4double finalE = kineticEnergy - gammaEnergy;
543
544  // stop tracking and create new secondary instead of primary
545  if(gammaEnergy > SecondaryThreshold()) {
546    fParticleChange->ProposeTrackStatus(fStopAndKill);
547    fParticleChange->SetProposedKineticEnergy(0.0);
548    G4DynamicParticle* el = 
549      new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
550                            direction, finalE);
551    vdp->push_back(el);
552
553    // continue tracking
554  } else {
555    fParticleChange->SetProposedMomentumDirection(direction);
556    fParticleChange->SetProposedKineticEnergy(finalE);
557  }
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
561
562
Note: See TracBrowser for help on using the repository browser.