Geant4  10.01.p01
G4eBremsstrahlungRelModel.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 // $Id: G4eBremsstrahlungRelModel.cc 83685 2014-09-09 12:39:00Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eBremsstrahlungRelModel
34 //
35 // Author: Andreas Schaelicke
36 //
37 // Creation date: 12.08.2008
38 //
39 // Modifications:
40 //
41 // 13.11.08 add SetLPMflag and SetLPMconstant methods
42 // 13.11.08 change default LPMconstant value
43 // 13.10.10 add angular distributon interface (VI)
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 
57 #include "G4PhysicalConstants.hh"
58 #include "G4SystemOfUnits.hh"
59 #include "G4Electron.hh"
60 #include "G4Positron.hh"
61 #include "G4Gamma.hh"
62 #include "Randomize.hh"
63 #include "G4Material.hh"
64 #include "G4Element.hh"
65 #include "G4ElementVector.hh"
66 #include "G4ProductionCutsTable.hh"
68 #include "G4LossTableManager.hh"
69 #include "G4ModifiedTsai.hh"
70 #include "G4DipBustGenerator.hh"
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73 
74 const G4double
75 G4eBremsstrahlungRelModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
76  0.5917, 0.7628, 0.8983, 0.9801 };
77 const G4double
78 G4eBremsstrahlungRelModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
79  0.1813, 0.1569, 0.1112, 0.0506 };
80 const G4double
81 G4eBremsstrahlungRelModel::Fel_light[] = {0., 5.31 , 4.79 , 4.74 , 4.71};
82 const G4double
83 G4eBremsstrahlungRelModel::Finel_light[] = {0., 6.144 , 5.621 , 5.805 , 5.924};
84 
85 using namespace std;
86 
88  const G4ParticleDefinition* p, const G4String& nam)
89  : G4VEmModel(nam),
90  particle(0),
91  bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
92  isElectron(true),
93  fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
94  fLPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)*0.5),
95  fXiLPM(0), fPhiLPM(0), fGLPM(0),
96  use_completescreening(false)
97 {
98  fParticleChange = 0;
100 
101  lowestKinEnergy = 1.0*MeV;
103 
105 
106  SetLPMFlag(true);
107  //SetAngularDistribution(new G4ModifiedTsai());
109 
112  = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
113 
114  energyThresholdLPM = 1.e39;
115 
117  if(p) { SetParticle(p); }
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
123 {
124  facFel = G4Log(184.15);
125  facFinel = G4Log(1194.);
126 
127  preS1 = 1./(184.15*184.15);
128  logTwo = G4Log(2.);
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132 
134 {
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138 
140 {
141  particle = p;
142  particleMass = p->GetPDGMass();
143  if(p == G4Electron::Electron()) { isElectron = true; }
144  else { isElectron = false;}
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
150  const G4Material* mat,
151  G4double kineticEnergy)
152 {
155 
156  // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
157  if (LPMFlag()) {
159  } else {
160  energyThresholdLPM=1.e39; // i.e. do not use LPM effect
161  }
162  // calculate threshold for density effect
163  kinEnergy = kineticEnergy;
164  totalEnergy = kineticEnergy + particleMass;
166 
167  // define critical gamma energies (important for integration/dicing)
168  klpm=totalEnergy*totalEnergy/lpmEnergy;
169  kp=sqrt(densityCorr);
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173 
175  const G4DataVector& cuts)
176 {
177  if(p) { SetParticle(p); }
178 
179  currentZ = 0.;
180 
181  if(IsMaster() && LowEnergyLimit() < HighEnergyLimit()) {
182  InitialiseElementSelectors(p, cuts);
183  }
184 
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189 
191  G4VEmModel* masterModel)
192 {
193  if(LowEnergyLimit() < HighEnergyLimit()) {
195  }
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199 
200 G4double
202  const G4ParticleDefinition*,
203  G4double cut)
204 {
205  return std::max(lowestKinEnergy, cut);
206 }
207 
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
209 
211  const G4Material* material,
212  const G4ParticleDefinition* p,
213  G4double kineticEnergy,
214  G4double cutEnergy)
215 {
216  if(!particle) { SetParticle(p); }
217  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
218  G4double cut = std::min(cutEnergy, kineticEnergy);
219  if(cut == 0.0) { return 0.0; }
220 
221  SetupForMaterial(particle, material,kineticEnergy);
222 
223  const G4ElementVector* theElementVector = material->GetElementVector();
224  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
225 
226  G4double dedx = 0.0;
227 
228  // loop for elements in the material
229  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
230 
231  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
232  SetCurrentElement((*theElementVector)[i]->GetZ());
233 
234  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
235  }
236  dedx *= bremFactor;
237 
238 
239  return dedx;
240 }
241 
242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
243 
245 {
246  G4double loss = 0.0;
247 
248  // number of intervals and integration step
249  G4double vcut = cut/totalEnergy;
250  G4int n = (G4int)(20*vcut) + 3;
251  G4double delta = vcut/G4double(n);
252 
253  G4double e0 = 0.0;
254  G4double xs;
255 
256  // integration
257  for(G4int l=0; l<n; l++) {
258 
259  for(G4int i=0; i<8; i++) {
260 
261  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
262 
265  } else {
266  xs = ComputeDXSectionPerAtom(eg);
267  }
268  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
269  }
270  e0 += delta;
271  }
272 
273  loss *= delta*totalEnergy;
274 
275  return loss;
276 }
277 
278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279 
281  const G4ParticleDefinition* p,
282  G4double kineticEnergy,
283  G4double Z, G4double,
284  G4double cutEnergy,
285  G4double maxEnergy)
286 {
287  if(!particle) { SetParticle(p); }
288  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
289  G4double cut = std::min(cutEnergy, kineticEnergy);
290  G4double tmax = std::min(maxEnergy, kineticEnergy);
291 
292  if(cut >= tmax) { return 0.0; }
293 
295 
296  G4double cross = ComputeXSectionPerAtom(cut);
297 
298  // allow partial integration
299  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
300 
301  cross *= Z*Z*bremFactor;
302 
303  return cross;
304 }
305 
306 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
307 
308 
310 {
311  G4double cross = 0.0;
312 
313  // number of intervals and integration step
314  G4double vcut = G4Log(cut/totalEnergy);
316  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
317  // n=1; // integration test
318  G4double delta = (vmax - vcut)/G4double(n);
319 
320  G4double e0 = vcut;
321  G4double xs;
322 
323  // integration
324  for(G4int l=0; l<n; l++) {
325 
326  for(G4int i=0; i<8; i++) {
327 
328  G4double eg = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
329 
332  } else {
333  xs = ComputeDXSectionPerAtom(eg);
334  }
335  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
336  }
337  e0 += delta;
338  }
339 
340  cross *= delta;
341 
342  return cross;
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
346 
348 {
349  // *** calculate lpm variable s & sprime ***
350  // Klein eqs. (78) & (79)
351  G4double sprime = sqrt(0.125*k*lpmEnergy/(totalEnergy*(totalEnergy-k)));
352 
353  G4double s1 = preS1*z23;
354  G4double logS1 = 2./3.*lnZ-2.*facFel;
355  G4double logTS1 = logTwo+logS1;
356 
357  xiLPM = 2.;
358 
359  if (sprime>1)
360  xiLPM = 1.;
361  else if (sprime>sqrt(2.)*s1) {
362  G4double h = G4Log(sprime)/logTS1;
363  xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
364  }
365 
366  G4double s0 = sprime/sqrt(xiLPM);
367 
368  // *** merging with density effect*** should be only necessary in region
369  // "close to" kp, e.g. k<100*kp using Ter-Mikaelian eq. (20.9)
370  G4double k2 = k*k;
371  s0 *= (1 + (densityCorr/k2) );
372 
373  // recalculate Xi using modified s above
374  // Klein eq. (75)
375  xiLPM = 1.;
376  if (s0<=s1) xiLPM = 2.;
377  else if ( (s1<s0) && (s0<=1) ) { xiLPM = 1. + G4Log(s0)/logS1; }
378 
379 
380  // *** calculate supression functions phi and G ***
381  // Klein eqs. (77)
382  G4double s2=s0*s0;
383  G4double s3=s0*s2;
384  G4double s4=s2*s2;
385 
386  if (s0<0.1) {
387  // high suppression limit
388  phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
389  - 57.69873135166053*s4;
390  gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
391  }
392  else if (s0<1.9516) {
393  // intermediate suppression
394  // using eq.77 approxim. valid s<2.
395  phiLPM = 1.-G4Exp(-6.*s0*(1.+(3.-pi)*s0)
396  +s3/(0.623+0.795*s0+0.658*s2));
397  if (s0<0.415827397755) {
398  // using eq.77 approxim. valid 0.07<s<2
399  G4double psiLPM = 1-G4Exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
400  gLPM = 3*psiLPM-2*phiLPM;
401  }
402  else {
403  // using alternative parametrisiation
404  G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
405  + s3*0.67282686077812381 + s4*-0.1207722909879257;
406  gLPM = tanh(pre);
407  }
408  }
409  else {
410  // low suppression limit valid s>2.
411  phiLPM = 1. - 0.0119048/s4;
412  gLPM = 1. - 0.0230655/s4;
413  }
414 
415  // *** make sure suppression is smaller than 1 ***
416  // *** caused by Migdal approximation in xi ***
417  if (xiLPM*phiLPM>1. || s0>0.57) { xiLPM=1./phiLPM; }
418 }
419 
420 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
421 
422 
424 // Ultra relativistic model
425 // only valid for very high energies, but includes LPM suppression
426 // * complete screening
427 {
428  if(gammaEnergy < 0.0) { return 0.0; }
429 
430  G4double y = gammaEnergy/totalEnergy;
431  G4double y2 = y*y*.25;
432  G4double yone2 = (1.-y+2.*y2);
433 
434  // ** form factors complete screening case **
435 
436  // ** calc LPM functions -- include ter-mikaelian merging with density effect **
437  // G4double xiLPM, gLPM, phiLPM; // to be made member variables !!!
438  CalcLPMFunctions(gammaEnergy);
439 
440  G4double mainLPM = xiLPM*(y2 * gLPM + yone2*phiLPM) * ( (Fel-fCoulomb) + Finel/currentZ );
441  G4double secondTerm = (1.-y)/12.*(1.+1./currentZ);
442 
443  G4double cross = mainLPM+secondTerm;
444  return cross;
445 }
446 
447 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
448 
450 // Relativistic model
451 // only valid for high energies (and if LPM suppression does not play a role)
452 // * screening according to thomas-fermi-Model (only valid for Z>5)
453 // * no LPM effect
454 {
455 
456  if(gammaEnergy < 0.0) { return 0.0; }
457 
458  G4double y = gammaEnergy/totalEnergy;
459 
460  G4double main=0.,secondTerm=0.;
461 
463  // ** form factors complete screening case **
464  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
465  secondTerm = (1.-y)/12.*(1.+1./currentZ);
466  }
467  else {
468  // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
469  G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
470  G4double gg=dd/z13;
471  G4double eps=dd/z23;
472  G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ);
473  G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ);
474 
475  main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
476  secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
477  }
478  G4double cross = main+secondTerm;
479  return cross;
480 }
481 
482 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
483 
485  std::vector<G4DynamicParticle*>* vdp,
486  const G4MaterialCutsCouple* couple,
487  const G4DynamicParticle* dp,
488  G4double cutEnergy,
489  G4double maxEnergy)
490 {
491  G4double kineticEnergy = dp->GetKineticEnergy();
492  if(kineticEnergy < LowEnergyLimit()) { return; }
493  G4double cut = std::min(cutEnergy, kineticEnergy);
494  G4double emax = std::min(maxEnergy, kineticEnergy);
495  if(cut >= emax) { return; }
496 
497  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
498 
499  const G4Element* elm =
500  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
501  SetCurrentElement(elm->GetZ());
502 
503  kinEnergy = kineticEnergy;
504  totalEnergy = kineticEnergy + particleMass;
506 
507  //G4double fmax= fMax;
508  G4bool highe = true;
509  if(totalEnergy < energyThresholdLPM) { highe = false; }
510 
511  G4double xmin = G4Log(cut*cut + densityCorr);
512  G4double xmax = G4Log(emax*emax + densityCorr);
513  G4double gammaEnergy, f, x;
514 
515  do {
516  x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
517  if(x < 0.0) { x = 0.0; }
518  gammaEnergy = sqrt(x);
519  if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
520  else { f = ComputeDXSectionPerAtom(gammaEnergy); }
521 
522  if ( f > fMax ) {
523  G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
524  << f << " > " << fMax
525  << " Egamma(MeV)= " << gammaEnergy
526  << " Ee(MeV)= " << kineticEnergy
527  << " " << GetName()
528  << G4endl;
529  }
530 
531  } while (f < fMax*G4UniformRand());
532 
533  //
534  // angles of the emitted gamma. ( Z - axis along the parent particle)
535  // use general interface
536  //
537 
538  G4ThreeVector gammaDirection =
539  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
540  G4lrint(currentZ),
541  couple->GetMaterial());
542 
543  // create G4DynamicParticle object for the Gamma
544  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
545  gammaEnergy);
546  vdp->push_back(gamma);
547 
548  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
549  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
550  - gammaEnergy*gammaDirection).unit();
551 
552  // energy of primary
553  G4double finalE = kineticEnergy - gammaEnergy;
554 
555  // stop tracking and create new secondary instead of primary
556  if(gammaEnergy > SecondaryThreshold()) {
559  G4DynamicParticle* el =
560  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
561  direction, finalE);
562  vdp->push_back(el);
563 
564  // continue tracking
565  } else {
568  }
569 }
570 
571 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
572 
573 
G4double Psi1(G4double, G4double)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:622
static const double MeV
Definition: G4SIunits.hh:193
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:650
G4bool isElectron(G4int ityp)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:615
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
const G4double pi
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:598
G4bool LPMFlag() const
Definition: G4VEmModel.hh:657
G4bool IsMaster() const
Definition: G4VEmModel.hh:699
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
static const G4double eps
G4double Psi1M2(G4double, G4double)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
static G4NistManager * Instance()
const G4ParticleDefinition * particle
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
static const G4double Finel_light[5]
G4double GetElectronDensity() const
Definition: G4Material.hh:217
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
G4double Phi1(G4double, G4double)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:783
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetRadlen() const
Definition: G4Material.hh:220
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
int main(int argc, char **argv)
Definition: genwindef.cc:359
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:791
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double Fel_light[5]
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:755
G4double Phi1M2(G4double, G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:605
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double ComputeXSectionPerAtom(G4double cutEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:776
G4ParticleChangeForLoss * fParticleChange
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetParticle(const G4ParticleDefinition *p)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:713
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:440
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:525
const G4Material * GetMaterial() const
G4double ComputeBremLoss(G4double cutEnergy)
void CalcLPMFunctions(G4double gammaEnergy)