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

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

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