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