Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 75582 2013-11-04 12:13:01Z 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),
92  isElectron(true),
95  fXiLPM(0), fPhiLPM(0), fGLPM(0),
96  use_completescreening(false),isInitialised(false)
97 {
98  fParticleChange = 0;
100 
101  lowestKinEnergy = 1.0*MeV;
102  SetLowEnergyLimit(lowestKinEnergy);
103 
105 
106  SetLPMFlag(true);
107  //SetAngularDistribution(new G4ModifiedTsai());
109 
110  particleMass = kinEnergy = totalEnergy = currentZ = z13 = z23 = lnZ = Fel
111  = Finel = fCoulomb = fMax = densityFactor = densityCorr = lpmEnergy
112  = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
113 
114  energyThresholdLPM = 1.e39;
115 
116  InitialiseConstants();
117  if(p) { SetParticle(p); }
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
122 void G4eBremsstrahlungRelModel::InitialiseConstants()
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 
139 void G4eBremsstrahlungRelModel::SetParticle(const G4ParticleDefinition* p)
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 {
153  densityFactor = mat->GetElectronDensity()*fMigdalConstant;
154  lpmEnergy = mat->GetRadlen()*fLPMconstant;
155 
156  // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
157  if (LPMFlag()) {
158  energyThresholdLPM=sqrt(densityFactor)*lpmEnergy;
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 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
176  const G4DataVector& cuts)
177 {
178  if(p) { SetParticle(p); }
179 
180  currentZ = 0.;
181 
182  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
183 
184  if(isInitialised) { return; }
186  isInitialised = true;
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190 
192  G4VEmModel* masterModel)
193 {
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198 
199 G4double
201  const G4ParticleDefinition*,
202  G4double cut)
203 {
204  return std::max(lowestKinEnergy, cut);
205 }
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208 
210  const G4Material* material,
211  const G4ParticleDefinition* p,
212  G4double kineticEnergy,
213  G4double cutEnergy)
214 {
215  if(!particle) { SetParticle(p); }
216  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
217  G4double cut = std::min(cutEnergy, kineticEnergy);
218  if(cut == 0.0) { return 0.0; }
219 
220  SetupForMaterial(particle, material,kineticEnergy);
221 
222  const G4ElementVector* theElementVector = material->GetElementVector();
223  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
224 
225  G4double dedx = 0.0;
226 
227  // loop for elements in the material
228  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
229 
230  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
231  SetCurrentElement((*theElementVector)[i]->GetZ());
232 
233  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
234  }
235  dedx *= bremFactor;
236 
237 
238  return dedx;
239 }
240 
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242 
243 G4double G4eBremsstrahlungRelModel::ComputeBremLoss(G4double cut)
244 {
245  G4double loss = 0.0;
246 
247  // number of intervals and integration step
248  G4double vcut = cut/totalEnergy;
249  G4int n = (G4int)(20*vcut) + 3;
250  G4double delta = vcut/G4double(n);
251 
252  G4double e0 = 0.0;
253  G4double xs;
254 
255  // integration
256  for(G4int l=0; l<n; l++) {
257 
258  for(G4int i=0; i<8; i++) {
259 
260  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
261 
262  if(totalEnergy > energyThresholdLPM) {
263  xs = ComputeRelDXSectionPerAtom(eg);
264  } else {
265  xs = ComputeDXSectionPerAtom(eg);
266  }
267  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
268  }
269  e0 += delta;
270  }
271 
272  loss *= delta*totalEnergy;
273 
274  return loss;
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
280  const G4ParticleDefinition* p,
281  G4double kineticEnergy,
282  G4double Z, G4double,
283  G4double cutEnergy,
284  G4double maxEnergy)
285 {
286  if(!particle) { SetParticle(p); }
287  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
288  G4double cut = std::min(cutEnergy, kineticEnergy);
289  G4double tmax = std::min(maxEnergy, kineticEnergy);
290 
291  if(cut >= tmax) { return 0.0; }
292 
294 
295  G4double cross = ComputeXSectionPerAtom(cut);
296 
297  // allow partial integration
298  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
299 
300  cross *= Z*Z*bremFactor;
301 
302  return cross;
303 }
304 
305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
306 
307 
308 G4double G4eBremsstrahlungRelModel::ComputeXSectionPerAtom(G4double cut)
309 {
310  G4double cross = 0.0;
311 
312  // number of intervals and integration step
313  G4double vcut = G4Log(cut/totalEnergy);
315  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
316  // n=1; // integration test
317  G4double delta = (vmax - vcut)/G4double(n);
318 
319  G4double e0 = vcut;
320  G4double xs;
321 
322  // integration
323  for(G4int l=0; l<n; l++) {
324 
325  for(G4int i=0; i<8; i++) {
326 
327  G4double eg = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
328 
329  if(totalEnergy > energyThresholdLPM) {
330  xs = ComputeRelDXSectionPerAtom(eg);
331  } else {
332  xs = ComputeDXSectionPerAtom(eg);
333  }
334  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
335  }
336  e0 += delta;
337  }
338 
339  cross *= delta;
340 
341  return cross;
342 }
343 
344 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
345 
346 void G4eBremsstrahlungRelModel::CalcLPMFunctions(G4double k)
347 {
348  // *** calculate lpm variable s & sprime ***
349  // Klein eqs. (78) & (79)
350  G4double sprime = sqrt(0.125*k*lpmEnergy/(totalEnergy*(totalEnergy-k)));
351 
352  G4double s1 = preS1*z23;
353  G4double logS1 = 2./3.*lnZ-2.*facFel;
354  G4double logTS1 = logTwo+logS1;
355 
356  xiLPM = 2.;
357 
358  if (sprime>1)
359  xiLPM = 1.;
360  else if (sprime>sqrt(2.)*s1) {
361  G4double h = G4Log(sprime)/logTS1;
362  xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
363  }
364 
365  G4double s0 = sprime/sqrt(xiLPM);
366 
367  // *** merging with density effect*** should be only necessary in region
368  // "close to" kp, e.g. k<100*kp using Ter-Mikaelian eq. (20.9)
369  G4double k2 = k*k;
370  s0 *= (1 + (densityCorr/k2) );
371 
372  // recalculate Xi using modified s above
373  // Klein eq. (75)
374  xiLPM = 1.;
375  if (s0<=s1) xiLPM = 2.;
376  else if ( (s1<s0) && (s0<=1) ) { xiLPM = 1. + G4Log(s0)/logS1; }
377 
378 
379  // *** calculate supression functions phi and G ***
380  // Klein eqs. (77)
381  G4double s2=s0*s0;
382  G4double s3=s0*s2;
383  G4double s4=s2*s2;
384 
385  if (s0<0.1) {
386  // high suppression limit
387  phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
388  - 57.69873135166053*s4;
389  gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
390  }
391  else if (s0<1.9516) {
392  // intermediate suppression
393  // using eq.77 approxim. valid s<2.
394  phiLPM = 1.-G4Exp(-6.*s0*(1.+(3.-pi)*s0)
395  +s3/(0.623+0.795*s0+0.658*s2));
396  if (s0<0.415827397755) {
397  // using eq.77 approxim. valid 0.07<s<2
398  G4double psiLPM = 1-G4Exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
399  gLPM = 3*psiLPM-2*phiLPM;
400  }
401  else {
402  // using alternative parametrisiation
403  G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
404  + s3*0.67282686077812381 + s4*-0.1207722909879257;
405  gLPM = tanh(pre);
406  }
407  }
408  else {
409  // low suppression limit valid s>2.
410  phiLPM = 1. - 0.0119048/s4;
411  gLPM = 1. - 0.0230655/s4;
412  }
413 
414  // *** make sure suppression is smaller than 1 ***
415  // *** caused by Migdal approximation in xi ***
416  if (xiLPM*phiLPM>1. || s0>0.57) { xiLPM=1./phiLPM; }
417 }
418 
419 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
420 
421 
422 G4double G4eBremsstrahlungRelModel::ComputeRelDXSectionPerAtom(G4double gammaEnergy)
423 // Ultra relativistic model
424 // only valid for very high energies, but includes LPM suppression
425 // * complete screening
426 {
427  if(gammaEnergy < 0.0) { return 0.0; }
428 
429  G4double y = gammaEnergy/totalEnergy;
430  G4double y2 = y*y*.25;
431  G4double yone2 = (1.-y+2.*y2);
432 
433  // ** form factors complete screening case **
434 
435  // ** calc LPM functions -- include ter-mikaelian merging with density effect **
436  // G4double xiLPM, gLPM, phiLPM; // to be made member variables !!!
437  CalcLPMFunctions(gammaEnergy);
438 
439  G4double mainLPM = xiLPM*(y2 * gLPM + yone2*phiLPM) * ( (Fel-fCoulomb) + Finel/currentZ );
440  G4double secondTerm = (1.-y)/12.*(1.+1./currentZ);
441 
442  G4double cross = mainLPM+secondTerm;
443  return cross;
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
449 // Relativistic model
450 // only valid for high energies (and if LPM suppression does not play a role)
451 // * screening according to thomas-fermi-Model (only valid for Z>5)
452 // * no LPM effect
453 {
454 
455  if(gammaEnergy < 0.0) { return 0.0; }
456 
457  G4double y = gammaEnergy/totalEnergy;
458 
459  G4double main=0.,secondTerm=0.;
460 
461  if (use_completescreening|| currentZ<5) {
462  // ** form factors complete screening case **
463  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
464  secondTerm = (1.-y)/12.*(1.+1./currentZ);
465  }
466  else {
467  // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
468  G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
469  G4double gg=dd/z13;
470  G4double eps=dd/z23;
471  G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ);
472  G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ);
473 
474  main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
475  secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
476  }
477  G4double cross = main+secondTerm;
478  return cross;
479 }
480 
481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
482 
484  std::vector<G4DynamicParticle*>* vdp,
485  const G4MaterialCutsCouple* couple,
486  const G4DynamicParticle* dp,
487  G4double cutEnergy,
488  G4double maxEnergy)
489 {
490  G4double kineticEnergy = dp->GetKineticEnergy();
491  if(kineticEnergy < LowEnergyLimit()) { return; }
492  G4double cut = std::min(cutEnergy, kineticEnergy);
493  G4double emax = std::min(maxEnergy, kineticEnergy);
494  if(cut >= emax) { return; }
495 
496  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
497 
498  const G4Element* elm =
499  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
500  SetCurrentElement(elm->GetZ());
501 
502  kinEnergy = kineticEnergy;
503  totalEnergy = kineticEnergy + particleMass;
505 
506  //G4double fmax= fMax;
507  G4bool highe = true;
508  if(totalEnergy < energyThresholdLPM) { highe = false; }
509 
510  G4double xmin = G4Log(cut*cut + densityCorr);
511  G4double xmax = G4Log(emax*emax + densityCorr);
512  G4double gammaEnergy, f, x;
513 
514  do {
515  x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
516  if(x < 0.0) { x = 0.0; }
517  gammaEnergy = sqrt(x);
518  if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
519  else { f = ComputeDXSectionPerAtom(gammaEnergy); }
520 
521  if ( f > fMax ) {
522  G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
523  << f << " > " << fMax
524  << " Egamma(MeV)= " << gammaEnergy
525  << " Ee(MeV)= " << kineticEnergy
526  << " " << GetName()
527  << G4endl;
528  }
529 
530  } while (f < fMax*G4UniformRand());
531 
532  //
533  // angles of the emitted gamma. ( Z - axis along the parent particle)
534  // use general interface
535  //
536 
537  G4ThreeVector gammaDirection =
538  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
539  G4lrint(currentZ),
540  couple->GetMaterial());
541 
542  // create G4DynamicParticle object for the Gamma
543  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
544  gammaEnergy);
545  vdp->push_back(gamma);
546 
547  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
548  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
549  - gammaEnergy*gammaDirection).unit();
550 
551  // energy of primary
552  G4double finalE = kineticEnergy - gammaEnergy;
553 
554  // stop tracking and create new secondary instead of primary
555  if(gammaEnergy > SecondaryThreshold()) {
558  G4DynamicParticle* el =
559  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
560  direction, finalE);
561  vdp->push_back(el);
562 
563  // continue tracking
564  } else {
567  }
568 }
569 
570 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
571 
572 
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:627
G4bool isElectron(G4int ityp)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
const char * p
Definition: xmltok.h:285
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4bool LPMFlag() const
Definition: G4VEmModel.hh:634
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
int main(int argc, char **argv)
Definition: genwindef.cpp:30
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
static G4NistManager * Instance()
const G4ParticleDefinition * particle
string material
Definition: eplot.py:19
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4double GetElectronDensity() const
Definition: G4Material.hh:215
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetRadlen() const
Definition: G4Material.hh:218
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
electron_Compton_length
Definition: hepunit.py:289
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:768
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
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:732
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:753
G4ParticleChangeForLoss * fParticleChange
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
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:690
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:433
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:510
const G4Material * GetMaterial() const