CMS 3D CMS Logo

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