Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eBremParametrizedModel.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: G4eBremParametrizedModel.cc 74581 2013-10-15 12:03:25Z gcosmo $
27 // GEANT4 tag $Name: geant4-09-04 $
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class file
32 //
33 //
34 // File name: G4eBremParametrizedModel
35 //
36 // Author: Andreas Schaelicke
37 //
38 // Creation date: 06.04.2011
39 //
40 // Modifications:
41 //
42 // Main References:
43 // - based on G4eBremsstrahlungModel and G4eBremsstrahlungRelModel
44 // -------------------------------------------------------------------
45 //
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "G4Electron.hh"
53 #include "G4Positron.hh"
54 #include "G4Gamma.hh"
55 #include "Randomize.hh"
56 #include "G4Material.hh"
57 #include "G4Element.hh"
58 #include "G4ElementVector.hh"
59 #include "G4ProductionCutsTable.hh"
61 #include "G4LossTableManager.hh"
62 #include "G4ModifiedTsai.hh"
63 #include "G4Exp.hh"
64 #include "G4Log.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
68 const G4double G4eBremParametrizedModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
69  0.5917, 0.7628, 0.8983, 0.9801 };
70 const G4double G4eBremParametrizedModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
71  0.1813, 0.1569, 0.1112, 0.0506 };
72 
73 static const G4double tlow = 1.*CLHEP::MeV;
74 
75 //
76 // GEANT4 internal units.
77 //
78 static const G4double
79  ah10 = 4.67733E+00, ah11 =-6.19012E-01, ah12 = 2.02225E-02,
80  ah20 =-7.34101E+00, ah21 = 1.00462E+00, ah22 =-3.20985E-02,
81  ah30 = 2.93119E+00, ah31 =-4.03761E-01, ah32 = 1.25153E-02;
82 
83 static const G4double
84  bh10 = 4.23071E+00, bh11 =-6.10995E-01, bh12 = 1.95531E-02,
85  bh20 =-7.12527E+00, bh21 = 9.69160E-01, bh22 =-2.74255E-02,
86  bh30 = 2.69925E+00, bh31 =-3.63283E-01, bh32 = 9.55316E-03;
87 
88 static const G4double
89  al00 =-2.05398E+00, al01 = 2.38815E-02, al02 = 5.25483E-04,
90  al10 =-7.69748E-02, al11 =-6.91499E-02, al12 = 2.22453E-03,
91  al20 = 4.06463E-02, al21 =-1.01281E-02, al22 = 3.40919E-04;
92 
93 static const G4double
94  bl00 = 1.04133E+00, bl01 =-9.43291E-03, bl02 =-4.54758E-04,
95  bl10 = 1.19253E-01, bl11 = 4.07467E-02, bl12 =-1.30718E-03,
96  bl20 =-1.59391E-02, bl21 = 7.27752E-03, bl22 =-1.94405E-04;
97 
98 using namespace std;
99 
101  const G4String& nam)
102  : G4VEmModel(nam),
103  particle(0),
104  isElectron(true),
107  isInitialised(false)
108 {
110 
111  minThreshold = 0.1*keV;
112  lowKinEnergy = 10.*MeV;
113  SetLowEnergyLimit(lowKinEnergy);
114 
116 
118 
120  = densityFactor = densityCorr = fMax = fCoulomb = 0.;
121 
122  InitialiseConstants();
123  if(p) { SetParticle(p); }
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
128 void G4eBremParametrizedModel::InitialiseConstants()
129 {
130  facFel = G4Log(184.15);
131  facFinel = G4Log(1194.);
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135 
137 {
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
142 void G4eBremParametrizedModel::SetParticle(const G4ParticleDefinition* p)
143 {
144  particle = p;
145  particleMass = p->GetPDGMass();
146  if(p == G4Electron::Electron()) { isElectron = true; }
147  else { isElectron = false;}
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
153  const G4MaterialCutsCouple*)
154 {
155  return minThreshold;
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
159 
161  const G4Material* mat,
162  G4double kineticEnergy)
163 {
164  densityFactor = mat->GetElectronDensity()*fMigdalConstant;
165 
166  // calculate threshold for density effect
167  kinEnergy = kineticEnergy;
168  totalEnergy = kineticEnergy + particleMass;
170 }
171 
172 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
176  const G4DataVector& cuts)
177 {
178  if(p) { SetParticle(p); }
179 
180  lowKinEnergy = LowEnergyLimit();
181 
182  currentZ = 0.;
183 
184  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
185 
186  if(isInitialised) { return; }
188  isInitialised = true;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
194  G4VEmModel* masterModel)
195 {
197 }
198 
199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200 
202  const G4Material* material,
203  const G4ParticleDefinition* p,
204  G4double kineticEnergy,
205  G4double cutEnergy)
206 {
207  if(!particle) { SetParticle(p); }
208  if(kineticEnergy < lowKinEnergy) { return 0.0; }
209  G4double cut = std::min(cutEnergy, kineticEnergy);
210  if(cut == 0.0) { return 0.0; }
211 
212  SetupForMaterial(particle, material,kineticEnergy);
213 
214  const G4ElementVector* theElementVector = material->GetElementVector();
215  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
216 
217  G4double dedx = 0.0;
218 
219  // loop for elements in the material
220  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
221 
222  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
223  SetCurrentElement((*theElementVector)[i]->GetZ());
224 
225  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
226  }
227  dedx *= bremFactor;
228 
229  return dedx;
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
233 
234 G4double G4eBremParametrizedModel::ComputeBremLoss(G4double cut)
235 {
236  G4double loss = 0.0;
237 
238  // number of intervals and integration step
239  G4double vcut = cut/totalEnergy;
240  G4int n = (G4int)(20*vcut) + 3;
241  G4double delta = vcut/G4double(n);
242 
243  G4double e0 = 0.0;
244  G4double xs;
245 
246  // integration
247  for(G4int l=0; l<n; l++) {
248 
249  for(G4int i=0; i<8; i++) {
250 
251  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
252 
253  xs = ComputeDXSectionPerAtom(eg);
254 
255  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
256  }
257  e0 += delta;
258  }
259 
260  loss *= delta*totalEnergy;
261 
262  return loss;
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266 
268  const G4ParticleDefinition* p,
269  G4double kineticEnergy,
270  G4double Z, G4double,
271  G4double cutEnergy,
272  G4double maxEnergy)
273 {
274  if(!particle) { SetParticle(p); }
275  if(kineticEnergy < lowKinEnergy) { return 0.0; }
276  G4double cut = std::min(cutEnergy, kineticEnergy);
277  G4double tmax = std::min(maxEnergy, kineticEnergy);
278 
279  if(cut >= tmax) { return 0.0; }
280 
281  SetCurrentElement(Z);
282 
283  G4double cross = ComputeXSectionPerAtom(cut);
284 
285  // allow partial integration
286  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
287 
288  cross *= Z*Z*bremFactor;
289 
290  return cross;
291 }
292 
293 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
294 
295 G4double G4eBremParametrizedModel::ComputeXSectionPerAtom(G4double cut)
296 {
297  G4double cross = 0.0;
298 
299  // number of intervals and integration step
300  G4double vcut = G4Log(cut/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 = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
315 
316  xs = ComputeDXSectionPerAtom(eg);
317 
318  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
319  }
320  e0 += delta;
321  }
322 
323  cross *= delta;
324 
325  return cross;
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329 
330 // compute the value of the screening function 3*PHI1 - PHI2
331 
332 G4double G4eBremParametrizedModel::ScreenFunction1(G4double ScreenVariable)
333 {
334  G4double screenVal;
335 
336  if (ScreenVariable > 1.)
337  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
338  else
339  screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
340 
341  return screenVal;
342 }
343 
344 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
345 
346 // compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
347 
348 G4double G4eBremParametrizedModel::ScreenFunction2(G4double ScreenVariable)
349 {
350  G4double screenVal;
351 
352  if (ScreenVariable > 1.)
353  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
354  else
355  screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
356 
357  return screenVal;
358 }
359 
360 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
361 
362 // Parametrized cross section
363 G4double G4eBremParametrizedModel::ComputeParametrizedDXSectionPerAtom(
364  G4double kineticEnergy,
365  G4double gammaEnergy, G4double Z)
366 {
367  SetCurrentElement(Z);
368  G4double FZ = lnZ* (4.- 0.55*lnZ);
369  G4double Z3 = z13;
370  G4double ZZ = z13*nist->GetZ13(G4lrint(Z)+1);
371 
372  totalEnergy = kineticEnergy + electron_mass_c2;
373 
374  // G4double x, epsil, greject, migdal, grejmax, q;
375  G4double epsil, greject;
376  G4double U = G4Log(kineticEnergy/electron_mass_c2);
377  G4double U2 = U*U;
378 
379  // precalculated parameters
380  G4double ah, bh;
381 
382  if (kineticEnergy > tlow) {
383 
384  G4double ah1 = ah10 + ZZ* (ah11 + ZZ* ah12);
385  G4double ah2 = ah20 + ZZ* (ah21 + ZZ* ah22);
386  G4double ah3 = ah30 + ZZ* (ah31 + ZZ* ah32);
387 
388  G4double bh1 = bh10 + ZZ* (bh11 + ZZ* bh12);
389  G4double bh2 = bh20 + ZZ* (bh21 + ZZ* bh22);
390  G4double bh3 = bh30 + ZZ* (bh31 + ZZ* bh32);
391 
392  ah = 1. + (ah1*U2 + ah2*U + ah3) / (U2*U);
393  bh = 0.75 + (bh1*U2 + bh2*U + bh3) / (U2*U);
394 
395  // limit of the screening variable
396  G4double screenfac =
397  136.*electron_mass_c2/(Z3*totalEnergy);
398 
399  epsil = gammaEnergy/totalEnergy; // epsil = x*kineticEnergy/totalEnergy;
400  G4double screenvar = screenfac*epsil/(1.0-epsil);
401  G4double F1 = max(ScreenFunction1(screenvar) - FZ ,0.);
402  G4double F2 = max(ScreenFunction2(screenvar) - FZ ,0.);
403 
404 
405  greject = (F1 - epsil* (ah*F1 - bh*epsil*F2))/8.; // 1./(42.392 - FZ);
406 
407  std::cout << " yy = "<<epsil<<std::endl;
408  std::cout << " F1/(...) "<<F1/(42.392 - FZ)<<std::endl;
409  std::cout << " F2/(...) "<<F2/(42.392 - FZ)<<std::endl;
410  std::cout << " (42.392 - FZ) " << (42.392 - FZ) <<std::endl;
411 
412  } else {
413 
414  G4double al0 = al00 + ZZ* (al01 + ZZ* al02);
415  G4double al1 = al10 + ZZ* (al11 + ZZ* al12);
416  G4double al2 = al20 + ZZ* (al21 + ZZ* al22);
417 
418  G4double bl0 = bl00 + ZZ* (bl01 + ZZ* bl02);
419  G4double bl1 = bl10 + ZZ* (bl11 + ZZ* bl12);
420  G4double bl2 = bl20 + ZZ* (bl21 + ZZ* bl22);
421 
422  ah = al0 + al1*U + al2*U2;
423  bh = bl0 + bl1*U + bl2*U2;
424 
425  G4double x=gammaEnergy/kineticEnergy;
426  greject=(1. + x* (ah + bh*x));
427 
428  /*
429  // Compute the maximum of the rejection function
430  grejmax = max(1. + xmin* (ah + bh*xmin), 1.+ah+bh);
431  G4double xm = -ah/(2.*bh);
432  if ( xmin < xm && xm < xmax) grejmax = max(grejmax, 1.+ xm* (ah + bh*xm));
433  */
434  }
435 
436  return greject;
437 }
438 
439 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
440 
441 G4double G4eBremParametrizedModel::ComputeDXSectionPerAtom(G4double gammaEnergy)
442 {
443 
444  if(gammaEnergy < 0.0) { return 0.0; }
445 
446  G4double y = gammaEnergy/totalEnergy;
447 
448  G4double main=0.;
449  //secondTerm=0.;
450 
451  // ** form factors complete screening case **
452  // only valid for high energies (and if LPM suppression does not play a role)
453  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
454  // secondTerm = (1.-y)/12.*(1.+1./currentZ);
455 
456  std::cout<<" F1(0) "<<ScreenFunction1(0.) <<std::endl;
457  std::cout<<" F1(0) "<<ScreenFunction2(0.) <<std::endl;
458  std::cout<<"Ekin = "<<kinEnergy<<std::endl;
459  std::cout<<"Z = "<<currentZ<<std::endl;
460  std::cout<<"main = "<<main<<std::endl;
461  std::cout<<" y = "<<y<<std::endl;
462  std::cout<<" Fel-fCoulomb "<< (Fel-fCoulomb) <<std::endl;
463 
464  G4double main2 = ComputeParametrizedDXSectionPerAtom(kinEnergy,gammaEnergy,currentZ);
465  std::cout<<"main2 = "<<main2<<std::endl;
466  std::cout<<"main2tot = "<<main2 * ( (Fel-fCoulomb) + Finel/currentZ )/(Fel-fCoulomb);
467 
468  G4double cross = main2; //main+secondTerm;
469  return cross;
470 }
471 
472 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473 
475  std::vector<G4DynamicParticle*>* vdp,
476  const G4MaterialCutsCouple* couple,
477  const G4DynamicParticle* dp,
478  G4double cutEnergy,
479  G4double maxEnergy)
480 {
481  G4double kineticEnergy = dp->GetKineticEnergy();
482  if(kineticEnergy < lowKinEnergy) { return; }
483  G4double cut = std::min(cutEnergy, kineticEnergy);
484  G4double emax = std::min(maxEnergy, kineticEnergy);
485  if(cut >= emax) { return; }
486 
487  SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
488 
489  const G4Element* elm =
490  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
491  SetCurrentElement(elm->GetZ());
492 
493  kinEnergy = kineticEnergy;
494  totalEnergy = kineticEnergy + particleMass;
496 
497  G4double xmin = G4Log(cut*cut + densityCorr);
498  G4double xmax = G4Log(emax*emax + densityCorr);
499  G4double gammaEnergy, f, x;
500 
501  do {
502  x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
503  if(x < 0.0) x = 0.0;
504  gammaEnergy = sqrt(x);
505  f = ComputeDXSectionPerAtom(gammaEnergy);
506 
507  if ( f > fMax ) {
508  G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
509  << f << " > " << fMax
510  << " Egamma(MeV)= " << gammaEnergy
511  << " E(mEV)= " << kineticEnergy
512  << G4endl;
513  }
514 
515  } while (f < fMax*G4UniformRand());
516 
517  //
518  // angles of the emitted gamma. ( Z - axis along the parent particle)
519  // use general interface
520  //
521  G4ThreeVector gammaDirection =
522  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
523  G4lrint(currentZ),
524  couple->GetMaterial());
525 
526  // create G4DynamicParticle object for the Gamma
527  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
528  gammaEnergy);
529  vdp->push_back(gamma);
530 
531  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
532  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
533  - gammaEnergy*gammaDirection).unit();
534 
535  // energy of primary
536  G4double finalE = kineticEnergy - gammaEnergy;
537 
538  // stop tracking and create new secondary instead of primary
539  if(gammaEnergy > SecondaryThreshold()) {
542  G4DynamicParticle* el =
543  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
544  direction, finalE);
545  vdp->push_back(el);
546 
547  // continue tracking
548  } else {
551  }
552 }
553 
554 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
555 
556 
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
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
G4ParticleChangeForLoss * fParticleChange
const char * p
Definition: xmltok.h:285
G4double GetZ13(G4double Z)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4eBremParametrizedModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremParam")
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
int main(int argc, char **argv)
Definition: genwindef.cpp:30
static const G4double xgi[8]
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
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
const G4ThreeVector & GetMomentumDirection() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
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)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
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
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
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
static const G4double wgi[8]
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:433
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
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)