Geant4  10.02.p01
G4DNAUeharaScreenedRutherfordElasticModel.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 
27 #include "G4PhysicalConstants.hh"
28 #include "G4SystemOfUnits.hh"
30 
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
32 
33 using namespace std;
34 
35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36 
38  const G4String& nam) :
39  G4VEmModel(nam), isInitialised(false)
40 {
41  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
42  fpWaterDensity = 0;
43 
44  killBelowEnergy = 9 * eV; // As in Geant4-DNA Screened Rutherford model
45  lowEnergyLimit = 0 * eV;
46  intermediateEnergyLimit = 200 * eV; // Switch between two final state models
47  highEnergyLimit = 10 * keV;
50 
51  verboseLevel = 0;
52  // Verbosity scale:
53  // 0 = nothing
54  // 1 = warning for energy non-conservation
55  // 2 = details of energy budget
56  // 3 = calculation of cross sections, file openings, sampling of atoms
57  // 4 = entering in methods
58 
59  if (verboseLevel > 0)
60  {
61  G4cout << "Screened Rutherford Elastic model is constructed " << G4endl<< "Energy range: "
62  << lowEnergyLimit / eV << " eV - "
63  << highEnergyLimit / MeV << " MeV"
64  << G4endl;
65  }
67 
68  // Selection of computation method
69  // We do not recommend "true" usage with the current cumul. proba. settings
70  fasterCode = false;
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
76 {
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82  const G4DataVector& /*cuts*/)
83 {
84 
85  if (verboseLevel > 3)
86  G4cout << "Calling G4DNAUeharaScreenedRutherfordElasticModel::Initialise()"
87  << G4endl;
88 
89  // Energy limits
90 
92  {
93  G4cout << "G4DNAUeharaScreenedRutherfordElasticModel: low energy limit increased from " <<
94  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
96  }
97 
99  {
100  G4cout << "G4DNAUeharaScreenedRutherfordElasticModel: high energy limit decreased from " <<
101  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
103  }
104 
105  // Constants for final state by Brenner & Zaider
106  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
107  // Added clear for MT
108 
109  betaCoeff.clear();
110  betaCoeff.push_back(7.51525);
111  betaCoeff.push_back(-0.41912);
112  betaCoeff.push_back(7.2017E-3);
113  betaCoeff.push_back(-4.646E-5);
114  betaCoeff.push_back(1.02897E-7);
115 
116  deltaCoeff.clear();
117  deltaCoeff.push_back(2.9612);
118  deltaCoeff.push_back(-0.26376);
119  deltaCoeff.push_back(4.307E-3);
120  deltaCoeff.push_back(-2.6895E-5);
121  deltaCoeff.push_back(5.83505E-8);
122 
123  gamma035_10Coeff.clear();
124  gamma035_10Coeff.push_back(-1.7013);
125  gamma035_10Coeff.push_back(-1.48284);
126  gamma035_10Coeff.push_back(0.6331);
127  gamma035_10Coeff.push_back(-0.10911);
128  gamma035_10Coeff.push_back(8.358E-3);
129  gamma035_10Coeff.push_back(-2.388E-4);
130 
131  gamma10_100Coeff.clear();
132  gamma10_100Coeff.push_back(-3.32517);
133  gamma10_100Coeff.push_back(0.10996);
134  gamma10_100Coeff.push_back(-4.5255E-3);
135  gamma10_100Coeff.push_back(5.8372E-5);
136  gamma10_100Coeff.push_back(-2.4659E-7);
137 
138  gamma100_200Coeff.clear();
139  gamma100_200Coeff.push_back(2.4775E-2);
140  gamma100_200Coeff.push_back(-2.96264E-5);
141  gamma100_200Coeff.push_back(-1.20655E-7);
142 
143  //
144 
145  if( verboseLevel>0 )
146  {
147  G4cout << "Screened Rutherford elastic model is initialized " << G4endl
148  << "Energy range: "
149  << LowEnergyLimit() / eV << " eV - "
150  << HighEnergyLimit() / MeV << " MeV"
151  << G4endl;
152  }
153 
154  // Initialize water density pointer
156 
157  if (isInitialised)
158  { return;}
160  isInitialised = true;
161 
162  }
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165 
167  const G4ParticleDefinition* particleDefinition,
168  G4double ekin,
169  G4double,
170  G4double)
171 {
172  if (verboseLevel > 3)
173  {
174  G4cout
175  << "Calling CrossSectionPerVolume() of G4DNAUeharaScreenedRutherfordElasticModel"
176  << G4endl;
177  }
178 
179  // Calculate total cross section for model
180 
181  G4double sigma=0;
182 
183  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
184 
185  if(waterDensity!= 0.0)
186  //if (material == nistwater || material->GetBaseMaterial() == nistwater)
187  {
188 
189  if (ekin < highEnergyLimit)
190  {
191  if (ekin < killBelowEnergy) return DBL_MAX;
192 
193  G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
194  G4double n = ScreeningFactor(ekin,z);
195  G4double crossSection = RutherfordCrossSection(ekin, z);
196  sigma = pi * crossSection / (n * (n + 1.));
197  }
198 
199  if (verboseLevel > 2)
200  {
201  G4cout << "__________________________________" << G4endl;
202  G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO START" << G4endl;
203  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
204  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
205  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
206  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
207  G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO END" << G4endl;
208  }
209 
210  }
211 
212  return sigma*material->GetAtomicNumDensityVector()[1];
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 
218  G4double z)
219 {
220  //
221  // e^4 / K + m_e c^2 \^2
222  // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
223  // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) /
224  //
225  // Where K is the electron non-relativistic kinetic energy
226  //
227  // NIM 155, pp. 145-156, 1978
228 
229  G4double length = (e_squared * (k + electron_mass_c2))
230  / (4 * pi * epsilon0 * k * (k + 2 * electron_mass_c2));
231  G4double cross = z * (z + 1) * length * length;
232 
233  return cross;
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
239  G4double z)
240 {
241  // From Phys Med Biol 37 (1992) 1841-1858
242  // Z=7.42 for water
243 
244  const G4double constK(1.7E-5);
245 
246  G4double beta2;
247  beta2 = 1. - 1. / ((1. + k / electron_mass_c2) * (1. + k / electron_mass_c2));
248 
249  G4double etaC;
250  if (k < 50 * keV)
251  etaC = 1.198;
252  else
253  etaC = 1.13 + 3.76 * (z * z / (137 * 137 * beta2));
254 
255  G4double numerator = etaC * constK * std::pow(z, 2. / 3.);
256 
257  k /= electron_mass_c2;
258 
259  G4double denominator = k * (2 + k);
260 
261  G4double value = 0.;
262  if (denominator > 0.)
263  value = numerator / denominator; // form 3
264 
265  return value;
266 }
267 
268 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
269 
271  G4DynamicParticle*>* /*fvect*/,
272  const G4MaterialCutsCouple* /*couple*/,
273  const G4DynamicParticle* aDynamicElectron,
274  G4double,
275  G4double)
276 {
277 
278  if (verboseLevel > 3)
279  {
280  G4cout
281  << "Calling SampleSecondaries() of G4DNAUeharaScreenedRutherfordElasticModel"
282  << G4endl;
283  }
284 
285  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
286 
287  if (electronEnergy0 < killBelowEnergy)
288  {
292  return;
293  }
294 
295  G4double cosTheta = 0.;
296 
297  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
298  {
299  if (electronEnergy0<intermediateEnergyLimit)
300  {
301  if (verboseLevel > 3) G4cout << "---> Using Brenner & Zaider model" << G4endl;
302  cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
303  }
304 
305  if (electronEnergy0>=intermediateEnergyLimit)
306  {
307  if (verboseLevel > 3) G4cout << "---> Using Screened Rutherford model" << G4endl;
308  G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
309  cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
310  }
311 
312  G4double phi = 2. * pi * G4UniformRand();
313 
314  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
315  G4ThreeVector xVers = zVers.orthogonal();
316  G4ThreeVector yVers = zVers.cross(xVers);
317 
318  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
319  G4double yDir = xDir;
320  xDir *= std::cos(phi);
321  yDir *= std::sin(phi);
322 
323  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
324 
326 
328  }
329 
330 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
333 
335 {
336  // d sigma_el 1 beta(K)
337  // ------------ (K) ~ --------------------------------- + ---------------------------------
338  // d Omega (1 + 2 gamma(K) - cos(theta))^2 (1 + 2 delta(K) + cos(theta))^2
339  //
340  // Maximum is < 1/(4 gamma(K)^2) + beta(K)/((2+2delta(K))^2)
341  //
342  // Phys. Med. Biol. 29 N.4 (1983) 443-447
343 
344  // gamma(K), beta(K) and delta(K) are polynomials with coefficients for energy measured in eV
345 
346  k /= eV;
347 
348  G4double beta = std::exp(CalculatePolynomial(k, betaCoeff));
349  G4double delta = std::exp(CalculatePolynomial(k, deltaCoeff));
350  G4double gamma;
351 
352  if (k > 100.)
353  {
355  // Only in this case it is not the exponent of the polynomial
356  } else
357  {
358  if (k > 10)
359  {
360  gamma = std::exp(CalculatePolynomial(k, gamma10_100Coeff));
361  } else
362  {
363  gamma = std::exp(CalculatePolynomial(k, gamma035_10Coeff));
364  }
365  }
366 
367  // ***** Original method
368 
369  if (!fasterCode)
370  {
371 
372  G4double oneOverMax = 1.
373  / (1. / (4. * gamma * gamma)
374  + beta / ((2. + 2. * delta) * (2. + 2. * delta)));
375 
376  G4double cosTheta = 0.;
377  G4double leftDenominator = 0.;
378  G4double rightDenominator = 0.;
379  G4double fCosTheta = 0.;
380 
381  do
382  {
383  cosTheta = 2. * G4UniformRand()- 1.;
384 
385  leftDenominator = (1. + 2.*gamma - cosTheta);
386  rightDenominator = (1. + 2.*delta + cosTheta);
387  if ( (leftDenominator * rightDenominator) != 0. )
388  {
389  fCosTheta = oneOverMax * (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator));
390  }
391  }
392  while (fCosTheta < G4UniformRand());
393 
394  return cosTheta;
395  }
396 
397  // ***** Alternative method using cumulative probability
398 
399  if (fasterCode)
400  {
401 
402  G4double cosTheta = -1;
403  G4double cumul = 0;
404  G4double value = 0;
405  G4double leftDenominator = 0.;
406  G4double rightDenominator = 0.;
407 
408  // Number of integration steps in the -1,1 range
409  G4int iMax=200;
410 
411  G4double random = G4UniformRand();
412 
413  // Cumulate differential cross section
414  for (G4int i=0; i<iMax; i++)
415  {
416  cosTheta = -1 + i*2./(iMax-1);
417  leftDenominator = (1. + 2.*gamma - cosTheta);
418  rightDenominator = (1. + 2.*delta + cosTheta);
419  if ( (leftDenominator * rightDenominator) != 0. )
420  {
421  cumul = cumul + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator));
422  }
423  }
424 
425  // Select cosTheta
426  for (G4int i=0; i<iMax; i++)
427  {
428  cosTheta = -1 + i*2./(iMax-1);
429  leftDenominator = (1. + 2.*gamma - cosTheta);
430  rightDenominator = (1. + 2.*delta + cosTheta);
431  if (cumul !=0 && (leftDenominator * rightDenominator) != 0.)
432  value = value + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator)) / cumul;
433  if (random < value) break;
434  }
435 
436  return cosTheta;
437  }
438 
439  return 0.;
440 
441 }
442 
443 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
444 
446  std::vector<
447  G4double>& vec)
448 {
449  // Sum_{i=0}^{size-1} vector_i k^i
450  //
451  // Phys. Med. Biol. 29 N.4 (1983) 443-447
452 
453  G4double result = 0.;
454  size_t size = vec.size();
455 
456  while (size > 0)
457  {
458  size--;
459 
460  result *= k;
461  result += vec[size];
462  }
463 
464  return result;
465 }
466 
467 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
468 
470  G4double z)
471 {
472 
473  // d sigma_el sigma_Ruth(K)
474  // ------------ (K) ~ -----------------------------
475  // d Omega (1 + 2 n(K) - cos(theta))^2
476  //
477  // We extract cos(theta) distributed as (1 + 2 n(K) - cos(theta))^-2
478  //
479  // 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)
480  //
481  // Phys. Med. Biol. 45 (2000) 3171-3194
482 
483  // ***** Original method
484 
485  if (!fasterCode)
486  {
487 
488  G4double n = ScreeningFactor(k, z);
489 
490  G4double oneOverMax = (4. * n * n);
491 
492  G4double cosTheta = 0.;
493  G4double fCosTheta;
494 
495  do
496  {
497  cosTheta = 2. * G4UniformRand()- 1.;
498  fCosTheta = (1 + 2.*n - cosTheta);
499  if (fCosTheta !=0.) fCosTheta = oneOverMax / (fCosTheta*fCosTheta);
500  }
501  while (fCosTheta < G4UniformRand());
502 
503  return cosTheta;
504  }
505 
506  // ***** Alternative method using cumulative probability
507  if (fasterCode)
508  {
509 
510  G4double cosTheta = -1;
511  G4double cumul = 0;
512  G4double value = 0;
513  G4double n = ScreeningFactor(k, z);
514  G4double fCosTheta;
515 
516  // Number of integration steps in the -1,1 range
517  G4int iMax=200;
518 
519  G4double random = G4UniformRand();
520 
521  // Cumulate differential cross section
522  for (G4int i=0; i<iMax; i++)
523  {
524  cosTheta = -1 + i*2./(iMax-1);
525  fCosTheta = (1 + 2.*n - cosTheta);
526  if (fCosTheta !=0.) cumul = cumul + 1./(fCosTheta*fCosTheta);
527  }
528 
529  // Select cosTheta
530  for (G4int i=0; i<iMax; i++)
531  {
532  cosTheta = -1 + i*2./(iMax-1);
533  fCosTheta = (1 + 2.*n - cosTheta);
534  if (cumul !=0.) value = value + (1./(fCosTheta*fCosTheta)) / cumul;
535  if (random < value) break;
536  }
537  return cosTheta;
538  }
539 
540  return 0.;
541 }
542 
static const double cm
Definition: G4SIunits.hh:118
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4double z
Definition: TRTMaterials.hh:39
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
size_t GetIndex() const
Definition: G4Material.hh:262
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4DNAUeharaScreenedRutherfordElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAUeharaScreenedRutherfordElasticModel")
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4int n
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
static const double pi
Definition: G4SIunits.hh:74
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4double CalculatePolynomial(G4double k, std::vector< G4double > &vec)
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
#define DBL_MAX
Definition: templates.hh:83
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134