Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAScreenedRutherfordElasticModel.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: G4DNAScreenedRutherfordElasticModel.cc 97520 2016-06-03 14:23:17Z gcosmo $
27 //
28 
30 #include "G4PhysicalConstants.hh"
31 #include "G4SystemOfUnits.hh"
33 #include "G4Exp.hh"
34 #include "G4Log.hh"
35 
36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37 
38 using namespace std;
39 //#define SR_VERBOSE
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42 
45  const G4String& nam) :
46  G4VEmModel(nam), isInitialised(false)
47 {
48  fpWaterDensity = 0;
49 
50  lowEnergyLimit = 0 * eV;
51  intermediateEnergyLimit = 200 * eV; // Switch between two final state models
52  highEnergyLimit = 1. * MeV;
53  SetLowEnergyLimit(lowEnergyLimit);
54  SetHighEnergyLimit(highEnergyLimit);
55 
56  verboseLevel = 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, sampling of atoms
62  // 4 = entering in methods
63 
64 #ifdef SR_VERBOSE
65  if (verboseLevel > 0)
66  {
67  G4cout << "Screened Rutherford Elastic model is constructed "
68  << G4endl
69  << "Energy range: "
70  << lowEnergyLimit / eV << " eV - "
71  << highEnergyLimit / MeV << " MeV"
72  << G4endl;
73  }
74 #endif
76 
77  // Selection of computation method
78  // We do not recommend "true" usage with the current cumul. proba. settings
79  fasterCode = false;
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85 {
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
92  const G4DataVector& /*cuts*/)
93 {
94 #ifdef SR_VERBOSE
95  if (verboseLevel > 3)
96  {
97  G4cout << "Calling G4DNAScreenedRutherfordElasticModel::Initialise()"
98  << G4endl;
99  }
100 #endif
101 
102  if(particle->GetParticleName() != "e-")
103  {
104  G4Exception ("*** WARNING: the G4DNAScreenedRutherfordElasticModel is not "
105  "intented to be used with another particle than the electron",
106  "",FatalException,"") ;
107  }
108 
109  // Energy limits
110  if (LowEnergyLimit() < 9*eV)
111  {
112  G4Exception("*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
113  "not validated below 9 eV",
114  "",JustWarning,"") ;
115  }
116 
117  if (HighEnergyLimit() > 1*MeV)
118  {
119  G4Exception("*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
120  "not validated above 1 MeV",
121  "",JustWarning,"") ;
122  }
123 
124  //
125 #ifdef SR_VERBOSE
126  if( verboseLevel>0 )
127  {
128  G4cout << "Screened Rutherford elastic model is initialized " << G4endl
129  << "Energy range: "
130  << LowEnergyLimit() / eV << " eV - "
131  << HighEnergyLimit() / MeV << " MeV"
132  << G4endl;
133  }
134 #endif
135 
136  if (isInitialised) { return; } // return here, prevent reinit consts + pointer
137 
138  // Initialize water density pointer
139  fpWaterDensity = G4DNAMolecularMaterial::Instance()->
140  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
141 
143  isInitialised = true;
144 
145  // Constants for final state by Brenner & Zaider
146  // note: if called after if(isInitialised) no need for clear and resetting
147  // the values at every call
148 
149  betaCoeff=
150  {
151  7.51525,
152  -0.41912,
153  7.2017E-3,
154  -4.646E-5,
155  1.02897E-7};
156 
157  deltaCoeff=
158  {
159  2.9612,
160  -0.26376,
161  4.307E-3,
162  -2.6895E-5,
163  5.83505E-8};
164 
165  gamma035_10Coeff =
166  {
167  -1.7013,
168  -1.48284,
169  0.6331,
170  -0.10911,
171  8.358E-3,
172  -2.388E-4};
173 
174  gamma10_100Coeff =
175  {
176  -3.32517,
177  0.10996,
178  -4.5255E-3,
179  5.8372E-5,
180  -2.4659E-7};
181 
182  gamma100_200Coeff =
183  {
184  2.4775E-2,
185  -2.96264E-5,
186  -1.20655E-7};
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190 
193 #ifdef SR_VERBOSE
194  const G4ParticleDefinition* particleDefinition,
195 #else
196  const G4ParticleDefinition*,
197 #endif
198  G4double ekin,
199  G4double,
200  G4double)
201 {
202 #ifdef SR_VERBOSE
203  if (verboseLevel > 3)
204  {
205  G4cout << "Calling CrossSectionPerVolume() of "
206  "G4DNAScreenedRutherfordElasticModel"
207  << G4endl;
208  }
209 #endif
210 
211  // Calculate total cross section for model
212 
213  G4double sigma=0.;
214  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
215 
216  if(waterDensity!= 0.0)
217  {
218  if(ekin < HighEnergyLimit() && ekin >= LowEnergyLimit())
219  {
220  G4double z = 10.;
221  G4double n = ScreeningFactor(ekin,z);
222  G4double crossSection = RutherfordCrossSection(ekin, z);
223  sigma = pi * crossSection / (n * (n + 1.));
224  }
225 
226 #ifdef SR_VERBOSE
227  if (verboseLevel > 2)
228  {
229  G4cout << "__________________________________" << G4endl;
230  G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO START"
231  << G4endl;
232  G4cout << "=== Kinetic energy(eV)=" << ekin/eV
233  << " particle : " << particleDefinition->GetParticleName()
234  << G4endl;
235  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
236  << G4endl;
237  G4cout << "=== Cross section per water molecule (cm^-1)="
238  << sigma*waterDensity/(1./cm) << G4endl;
239  G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO END"
240  << G4endl;
241  }
242 #endif
243  }
244 
245  return sigma*waterDensity;
246 }
247 
248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
249 
250 G4double G4DNAScreenedRutherfordElasticModel::RutherfordCrossSection(G4double k,
251  G4double z)
252 {
253  //
254  // e^4 / K + m_e c^2 \^2
255  // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
256  // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) /
257  //
258  // Where K is the electron non-relativistic kinetic energy
259  //
260  // NIM 155, pp. 145-156, 1978
261 
262  G4double length = (e_squared * (k + electron_mass_c2))
263  / (4 * pi * epsilon0 * k * (k + 2 * electron_mass_c2));
264  G4double cross = z * (z + 1) * length * length;
265 
266  return cross;
267 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270 
271 G4double G4DNAScreenedRutherfordElasticModel::ScreeningFactor(G4double k,
272  G4double z)
273 {
274  //
275  // alpha_1 + beta_1 ln(K/eV) constK Z^(2/3)
276  // n(T) = -------------------------- -----------------
277  // K/(m_e c^2) 2 + K/(m_e c^2)
278  //
279  // Where K is the electron non-relativistic kinetic energy
280  //
281  // n(T) > 0 for T < ~ 400 MeV
282  //
283  // NIM 155, pp. 145-156, 1978
284  // Formulae (2) and (5)
285 
286  const G4double alpha_1(1.64);
287  const G4double beta_1(-0.0825);
288  const G4double constK(1.7E-5);
289 
290  G4double numerator = (alpha_1 + beta_1 * G4Log(k / eV)) * constK
291  * std::pow(z, 2. / 3.);
292 
293  k /= electron_mass_c2;
294 
295  G4double denominator = k * (2 + k);
296 
297  G4double value = 0.;
298  if (denominator > 0.) value = numerator / denominator;
299 
300  return value;
301 }
302 
303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
304 
306 SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
307  const G4MaterialCutsCouple* /*couple*/,
308  const G4DynamicParticle* aDynamicElectron,
309  G4double,
310  G4double)
311 {
312 #ifdef SR_VERBOSE
313  if (verboseLevel > 3)
314  {
315  G4cout << "Calling SampleSecondaries() of "
316  "G4DNAScreenedRutherfordElasticModel"
317  << G4endl;
318  }
319 #endif
320 
321  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
322  G4double cosTheta = 0.;
323 
324  // if (electronEnergy0 < highEnergyLimit)
325  {
326  if (electronEnergy0<intermediateEnergyLimit)
327  {
328 #ifdef SR_VERBOSE
329  if (verboseLevel > 3)
330  {G4cout << "---> Using Brenner & Zaider model" << G4endl;}
331 #endif
332  cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
333  }
334 
335  if (electronEnergy0>=intermediateEnergyLimit)
336  {
337 #ifdef SR_VERBOSE
338  if (verboseLevel > 3)
339  {G4cout << "---> Using Screened Rutherford model" << G4endl;}
340 #endif
341  G4double z = 10.;
342  cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
343  }
344 
345  G4double phi = 2. * pi * G4UniformRand();
346 
347  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
348  G4ThreeVector xVers = zVers.orthogonal();
349  G4ThreeVector yVers = zVers.cross(xVers);
350 
351  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
352  G4double yDir = xDir;
353  xDir *= std::cos(phi);
354  yDir *= std::sin(phi);
355 
356  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
357 
359 
361  }
362 }
363 
364 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
365 
366 G4double G4DNAScreenedRutherfordElasticModel::
367 BrennerZaiderRandomizeCosTheta(G4double k)
368 {
369  // d sigma_el 1 beta(K)
370  // ------------ (K) ~ --------------------------------- + ---------------------------------
371  // d Omega (1 + 2 gamma(K) - cos(theta))^2 (1 + 2 delta(K) + cos(theta))^2
372  //
373  // Maximum is < 1/(4 gamma(K)^2) + beta(K)/((2+2delta(K))^2)
374  //
375  // Phys. Med. Biol. 29 N.4 (1983) 443-447
376 
377  // gamma(K), beta(K) and delta(K) are polynomials with coefficients for energy measured in eV
378 
379  k /= eV;
380 
381  G4double beta = G4Exp(CalculatePolynomial(k, betaCoeff));
382  G4double delta = G4Exp(CalculatePolynomial(k, deltaCoeff));
383  G4double gamma;
384 
385  if (k > 100.)
386  {
387  gamma = CalculatePolynomial(k, gamma100_200Coeff);
388  // Only in this case it is not the exponent of the polynomial
389  }
390  else
391  {
392  if (k > 10)
393  {
394  gamma = G4Exp(CalculatePolynomial(k, gamma10_100Coeff));
395  }
396  else
397  {
398  gamma = G4Exp(CalculatePolynomial(k, gamma035_10Coeff));
399  }
400  }
401 
402  // ***** Original method
403 
404  if (!fasterCode)
405  {
406 
407  G4double oneOverMax = 1.
408  / (1. / (4. * gamma * gamma) + beta
409  / ((2. + 2. * delta) * (2. + 2. * delta)));
410 
411  G4double cosTheta = 0.;
412  G4double leftDenominator = 0.;
413  G4double rightDenominator = 0.;
414  G4double fCosTheta = 0.;
415 
416  do
417  {
418  cosTheta = 2. * G4UniformRand()- 1.;
419 
420  leftDenominator = (1. + 2.*gamma - cosTheta);
421  rightDenominator = (1. + 2.*delta + cosTheta);
422  if ( (leftDenominator * rightDenominator) != 0. )
423  {
424  fCosTheta = oneOverMax * (1./(leftDenominator*leftDenominator)
425  + beta/(rightDenominator*rightDenominator));
426  }
427  }
428  while (fCosTheta < G4UniformRand());
429 
430  return cosTheta;
431  }
432 
433  // ***** Alternative method using cumulative probability
434 
435  if (fasterCode)
436  {
437 
438  //
439  // modified by Shogo OKADA @ KEK, JP, 2016.2.27(Sat.)
440  //
441  // An integral of differential cross-section formula shown above this member function
442  // (integral variable: cos(theta), integral interval: [-1, x]) is as follows:
443  //
444  // 1.0 + x beta * (1 + x)
445  // I = --------------------- + ---------------------- (1)
446  // (a - x) * (a + 1.0) (b + x) * (b - 1.0)
447  //
448  // where a = 1.0 + 2.0 * gamma(K), b = 1.0 + 2.0 * delta(K)
449  //
450  // Then, a cumulative probability (cp) is as follows:
451  //
452  // cp 1.0 + x beta * (1 + x)
453  // ---- = --------------------- + ---------------------- (2)
454  // S (a - x) * (a + 1.0) (b + x) * (b - 1.0)
455  //
456  // where 1/S is the integral of differnetical cross-section (1) on interval [-1, 1]
457  //
458  // 1 2.0 2.0 * beta
459  // --- = ----------------------- + ----------------------- (3)
460  // S (a - 1.0) * (a + 1.0) (b + 1.0) * (b - 1.0)
461  //
462  // x is calculated from the quadratic equation derived from (2) and (3):
463  //
464  // A * x^2 + B * x + C = 0
465  //
466  // where A, B, anc C are coefficients of the equation:
467  // A = S * {(b - 1.0) - beta * (a + 1.0)} + cp * (a + 1.0) * (b - 1.0),
468  // B = S * {(b - 1.0) * (b + 1.0) + beta * (a - 1.0) * (a + 1.0)} - cp * (a + 1.0) * (b - 1.0) * (a - b)
469  // C = S * {b * (b - 1.0) + beta * a * (a + 1.0)} - cp * (a + 1.0) * (b - 1.0) * ab
470  //
471 
472  // sampling cumulative probability
474 
475  G4double a = 1.0 + 2.0 * gamma;
476  G4double b = 1.0 + 2.0 * delta;
477  G4double a1 = a - 1.0;
478  G4double a2 = a + 1.0;
479  G4double b1 = b - 1.0;
480  G4double b2 = b + 1.0;
481  G4double c1 = a - b;
482  G4double c2 = a * b;
483 
484  G4double S = 2.0 / (a1 * a2) + 2.0 * beta / (b1 * b2); S = 1.0 / S;
485 
486  // coefficients for the quadratic equation
487  G4double A = S * (b1 - beta * a2) + cp * a2 * b1;
488  G4double B = S * (b1 * b2 + beta * a1 * a2) - cp * a2 * b1 * c1;
489  G4double C = S * (b * b1 + beta * a * a2) - cp * a2 * b1 * c2;
490 
491  // calculate cos(theta)
492  return (-1.0 * B + std::sqrt(B * B - 4.0 * A * C)) / (2.0 * A);
493 
494  /*
495  G4double cosTheta = -1;
496  G4double cumul = 0;
497  G4double value = 0;
498  G4double leftDenominator = 0.;
499  G4double rightDenominator = 0.;
500 
501  // Number of integration steps in the -1,1 range
502  G4int iMax=200;
503 
504  G4double random = G4UniformRand();
505 
506  // Cumulate differential cross section
507  for (G4int i=0; i<iMax; i++)
508  {
509  cosTheta = -1 + i*2./(iMax-1);
510  leftDenominator = (1. + 2.*gamma - cosTheta);
511  rightDenominator = (1. + 2.*delta + cosTheta);
512  if ( (leftDenominator * rightDenominator) != 0. )
513  {
514  cumul = cumul + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator));
515  }
516  }
517 
518  // Select cosTheta
519  for (G4int i=0; i<iMax; i++)
520  {
521  cosTheta = -1 + i*2./(iMax-1);
522  leftDenominator = (1. + 2.*gamma - cosTheta);
523  rightDenominator = (1. + 2.*delta + cosTheta);
524  if (cumul !=0 && (leftDenominator * rightDenominator) != 0.)
525  value = value + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator)) / cumul;
526  if (random < value) break;
527  }
528 
529  return cosTheta;
530  */
531  }
532 
533  return 0.;
534 }
535 
536 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
537 
538 G4double G4DNAScreenedRutherfordElasticModel::
539 CalculatePolynomial(G4double k,
540  std::vector<G4double>& vec)
541 {
542  // Sum_{i=0}^{size-1} vector_i k^i
543  //
544  // Phys. Med. Biol. 29 N.4 (1983) 443-447
545 
546  G4double result = 0.;
547  size_t size = vec.size();
548 
549  while (size > 0)
550  {
551  size--;
552 
553  result *= k;
554  result += vec[size];
555  }
556 
557  return result;
558 }
559 
560 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
561 
562 G4double G4DNAScreenedRutherfordElasticModel::
563 ScreenedRutherfordRandomizeCosTheta(G4double k,
564  G4double z)
565 {
566 
567  // d sigma_el sigma_Ruth(K)
568  // ------------ (K) ~ -----------------------------
569  // d Omega (1 + 2 n(K) - cos(theta))^2
570  //
571  // We extract cos(theta) distributed as (1 + 2 n(K) - cos(theta))^-2
572  //
573  // Maximum is for theta=0: 1/(4 n(K)^2) (When n(K) is positive, that is always satisfied within the validity of the process)
574  //
575  // Phys. Med. Biol. 45 (2000) 3171-3194
576 
577  // ***** Original method
578 
579  if (!fasterCode)
580  {
581 
582  G4double n = ScreeningFactor(k, z);
583 
584  G4double oneOverMax = (4. * n * n);
585 
586  G4double cosTheta = 0.;
587  G4double fCosTheta;
588 
589  do
590  {
591  cosTheta = 2. * G4UniformRand()- 1.;
592  fCosTheta = (1 + 2.*n - cosTheta);
593  if (fCosTheta !=0.) fCosTheta = oneOverMax / (fCosTheta*fCosTheta);
594  }
595  while (fCosTheta < G4UniformRand());
596 
597  return cosTheta;
598  }
599 
600  // ***** Alternative method using cumulative probability
601  if (fasterCode)
602  {
603 
604  //
605  // modified by Shogo OKADA @ KEK, JP, 2016.2.27(Sat.)
606  //
607  // The cumulative probability (cp) is calculated by integrating
608  // the differential cross-section fomula with cos(theta):
609  //
610  // n(K) * (1.0 + cos(theta))
611  // cp = ---------------------------------
612  // 1.0 + 2.0 * n(K) - cos(theta)
613  //
614  // Then, cos(theta) is as follows:
615  //
616  // cp * (1.0 + 2.0 * n(K)) - n(K)
617  // cos(theta) = --------------------------------
618  // n(k) + cp
619  //
620  // where, K is kinetic energy, n(K) is screeing factor, and cp is cumulative probability
621  //
622 
623  G4double n = ScreeningFactor(k, z);
624  G4double cp = G4UniformRand();
625  G4double numerator = cp * (1.0 + 2.0 * n) - n;
626  G4double denominator = n + cp;
627  return numerator / denominator;
628 
629  /*
630  G4double cosTheta = -1;
631  G4double cumul = 0;
632  G4double value = 0;
633  G4double n = ScreeningFactor(k, z);
634  G4double fCosTheta;
635 
636  // Number of integration steps in the -1,1 range
637  G4int iMax=200;
638 
639  G4double random = G4UniformRand();
640 
641  // Cumulate differential cross section
642  for (G4int i=0; i<iMax; i++)
643  {
644  cosTheta = -1 + i*2./(iMax-1);
645  fCosTheta = (1 + 2.*n - cosTheta);
646  if (fCosTheta !=0.) cumul = cumul + 1./(fCosTheta*fCosTheta);
647  }
648 
649  // Select cosTheta
650  for (G4int i=0; i<iMax; i++)
651  {
652  cosTheta = -1 + i*2./(iMax-1);
653  fCosTheta = (1 + 2.*n - cosTheta);
654  if (cumul !=0.) value = value + (1./(fCosTheta*fCosTheta)) / cumul;
655  if (random < value) break;
656  }
657  return cosTheta;
658  */
659  }
660 
661  return 0.;
662 }
663 
664 
G4double G4ParticleHPJENDLHEData::G4double result
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
double S(double temp)
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
static constexpr double epsilon0
size_t GetIndex() const
Definition: G4Material.hh:262
double C(double temp)
double B(double temperature)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
static constexpr double electron_mass_c2
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4ThreeVector & GetMomentumDirection() const
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double e_squared
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4DNAMolecularMaterial * Instance()
Hep3Vector unit() const
Hep3Vector orthogonal() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4DNAScreenedRutherfordElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAScreenedRutherfordElasticModel")
static constexpr double pi
Definition: G4SIunits.hh:75
Hep3Vector cross(const Hep3Vector &) const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static const G4double cp
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132