Geant4  10.02
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 92074 2015-08-17 07:03:46Z 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  // Selection of computation method
72  // We do not recommend "true" usage with the current cumul. proba. settings
73  fasterCode = false;
74 }
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
79 {
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85  const G4DataVector& /*cuts*/)
86 {
87 
88  if (verboseLevel > 3)
89  G4cout << "Calling G4DNAScreenedRutherfordElasticModel::Initialise()"
90  << G4endl;
91 
92  // Energy limits
93 
95  {
96  G4cout << "G4DNAScreenedRutherfordElasticModel: low energy limit increased from " <<
97  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
99  }
100 
102  {
103  G4cout << "G4DNAScreenedRutherfordElasticModel: high energy limit decreased from " <<
104  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
106  }
107 
108  // Constants for final state by Brenner & Zaider
109  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
110  // Added clear for MT
111 
112  betaCoeff.clear();
113  betaCoeff.push_back(7.51525);
114  betaCoeff.push_back(-0.41912);
115  betaCoeff.push_back(7.2017E-3);
116  betaCoeff.push_back(-4.646E-5);
117  betaCoeff.push_back(1.02897E-7);
118 
119  deltaCoeff.clear();
120  deltaCoeff.push_back(2.9612);
121  deltaCoeff.push_back(-0.26376);
122  deltaCoeff.push_back(4.307E-3);
123  deltaCoeff.push_back(-2.6895E-5);
124  deltaCoeff.push_back(5.83505E-8);
125 
126  gamma035_10Coeff.clear();
127  gamma035_10Coeff.push_back(-1.7013);
128  gamma035_10Coeff.push_back(-1.48284);
129  gamma035_10Coeff.push_back(0.6331);
130  gamma035_10Coeff.push_back(-0.10911);
131  gamma035_10Coeff.push_back(8.358E-3);
132  gamma035_10Coeff.push_back(-2.388E-4);
133 
134  gamma10_100Coeff.clear();
135  gamma10_100Coeff.push_back(-3.32517);
136  gamma10_100Coeff.push_back(0.10996);
137  gamma10_100Coeff.push_back(-4.5255E-3);
138  gamma10_100Coeff.push_back(5.8372E-5);
139  gamma10_100Coeff.push_back(-2.4659E-7);
140 
141  gamma100_200Coeff.clear();
142  gamma100_200Coeff.push_back(2.4775E-2);
143  gamma100_200Coeff.push_back(-2.96264E-5);
144  gamma100_200Coeff.push_back(-1.20655E-7);
145 
146  //
147 
148  if( verboseLevel>0 )
149  {
150  G4cout << "Screened Rutherford elastic model is initialized " << G4endl
151  << "Energy range: "
152  << LowEnergyLimit() / eV << " eV - "
153  << HighEnergyLimit() / MeV << " MeV"
154  << G4endl;
155  }
156 
157  // Initialize water density pointer
159 
160  if (isInitialised)
161  { return;}
163  isInitialised = true;
164 
165 }
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168 
170  const G4ParticleDefinition* particleDefinition,
171  G4double ekin,
172  G4double,
173  G4double)
174 {
175  if (verboseLevel > 3)
176  G4cout << "Calling CrossSectionPerVolume() of G4DNAScreenedRutherfordElasticModel"
177  << G4endl;
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 
192  if (ekin < killBelowEnergy) return DBL_MAX;
193 
194  G4double z = 10.;
195  G4double n = ScreeningFactor(ekin,z);
196  G4double crossSection = RutherfordCrossSection(ekin, z);
197  sigma = pi * crossSection / (n * (n + 1.));
198  }
199 
200  if (verboseLevel > 2)
201  {
202  G4cout << "__________________________________" << G4endl;
203  G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO START" << G4endl;
204  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
205  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
206  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
207  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
208  G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO END" << G4endl;
209  }
210 
211  }
212 
213  return sigma*waterDensity;
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219  G4double z)
220 {
221  //
222  // e^4 / K + m_e c^2 \^2
223  // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
224  // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) /
225  //
226  // Where K is the electron non-relativistic kinetic energy
227  //
228  // NIM 155, pp. 145-156, 1978
229 
230  G4double length = (e_squared * (k + electron_mass_c2))
231  / (4 * pi * epsilon0 * k * (k + 2 * electron_mass_c2));
232  G4double cross = z * (z + 1) * length * length;
233 
234  return cross;
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
238 
240  G4double z)
241 {
242  //
243  // alpha_1 + beta_1 ln(K/eV) constK Z^(2/3)
244  // n(T) = -------------------------- -----------------
245  // K/(m_e c^2) 2 + K/(m_e c^2)
246  //
247  // Where K is the electron non-relativistic kinetic energy
248  //
249  // n(T) > 0 for T < ~ 400 MeV
250  //
251  // NIM 155, pp. 145-156, 1978
252  // Formulae (2) and (5)
253 
254  const G4double alpha_1(1.64);
255  const G4double beta_1(-0.0825);
256  const G4double constK(1.7E-5);
257 
258  G4double numerator = (alpha_1 + beta_1 * std::log(k / eV)) * constK
259  * std::pow(z, 2. / 3.);
260 
261  k /= electron_mass_c2;
262 
263  G4double denominator = k * (2 + k);
264 
265  G4double value = 0.;
266  if (denominator > 0.) value = numerator / denominator;
267 
268  return value;
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272 
274  G4DynamicParticle*>* /*fvect*/,
275  const G4MaterialCutsCouple* /*couple*/,
276  const G4DynamicParticle* aDynamicElectron,
277  G4double,
278  G4double)
279 {
280 
281  if (verboseLevel > 3)
282  {
283  G4cout << "Calling SampleSecondaries() of G4DNAScreenedRutherfordElasticModel"
284  << G4endl;
285  }
286 
287  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
288 
289  if (electronEnergy0 < killBelowEnergy)
290  {
294  return;
295  }
296 
297  G4double cosTheta = 0.;
298 
299  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
300  {
301  if (electronEnergy0<intermediateEnergyLimit)
302  {
303  if (verboseLevel > 3) G4cout << "---> Using Brenner & Zaider model" << G4endl;
304  cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
305  }
306 
307  if (electronEnergy0>=intermediateEnergyLimit)
308  {
309  if (verboseLevel > 3) G4cout << "---> Using Screened Rutherford model" << G4endl;
310  G4double z = 10.;
311  cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
312  }
313 
314  G4double phi = 2. * pi * G4UniformRand();
315 
316  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
317  G4ThreeVector xVers = zVers.orthogonal();
318  G4ThreeVector yVers = zVers.cross(xVers);
319 
320  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
321  G4double yDir = xDir;
322  xDir *= std::cos(phi);
323  yDir *= std::sin(phi);
324 
325  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
326 
328 
330  }
331 
332 }
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
335 
337 {
338  // d sigma_el 1 beta(K)
339  // ------------ (K) ~ --------------------------------- + ---------------------------------
340  // d Omega (1 + 2 gamma(K) - cos(theta))^2 (1 + 2 delta(K) + cos(theta))^2
341  //
342  // Maximum is < 1/(4 gamma(K)^2) + beta(K)/((2+2delta(K))^2)
343  //
344  // Phys. Med. Biol. 29 N.4 (1983) 443-447
345 
346  // gamma(K), beta(K) and delta(K) are polynomials with coefficients for energy measured in eV
347 
348  k /= eV;
349 
350  G4double beta = std::exp(CalculatePolynomial(k, betaCoeff));
351  G4double delta = std::exp(CalculatePolynomial(k, deltaCoeff));
352  G4double gamma;
353 
354  if (k > 100.)
355  {
357  // Only in this case it is not the exponent of the polynomial
358  }
359  else
360  {
361  if (k > 10)
362  {
363  gamma = std::exp(CalculatePolynomial(k, gamma10_100Coeff));
364  }
365  else
366  {
367  gamma = std::exp(CalculatePolynomial(k, gamma035_10Coeff));
368  }
369  }
370 
371  // ***** Original method
372 
373  if (!fasterCode)
374  {
375 
376  G4double oneOverMax = 1.
377  / (1. / (4. * gamma * gamma) + beta
378  / ((2. + 2. * delta) * (2. + 2. * delta)));
379 
380  G4double cosTheta = 0.;
381  G4double leftDenominator = 0.;
382  G4double rightDenominator = 0.;
383  G4double fCosTheta = 0.;
384 
385  do
386  {
387  cosTheta = 2. * G4UniformRand()- 1.;
388 
389  leftDenominator = (1. + 2.*gamma - cosTheta);
390  rightDenominator = (1. + 2.*delta + cosTheta);
391  if ( (leftDenominator * rightDenominator) != 0. )
392  {
393  fCosTheta = oneOverMax * (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator));
394  }
395  }
396  while (fCosTheta < G4UniformRand());
397 
398  return cosTheta;
399  }
400 
401  // ***** Alternative method using cumulative probability
402 
403  if (fasterCode)
404  {
405 
406  G4double cosTheta = -1;
407  G4double cumul = 0;
408  G4double value = 0;
409  G4double leftDenominator = 0.;
410  G4double rightDenominator = 0.;
411 
412  // Number of integration steps in the -1,1 range
413  G4int iMax=200;
414 
415  G4double random = G4UniformRand();
416 
417  // Cumulate differential cross section
418  for (G4int i=0; i<iMax; i++)
419  {
420  cosTheta = -1 + i*2./(iMax-1);
421  leftDenominator = (1. + 2.*gamma - cosTheta);
422  rightDenominator = (1. + 2.*delta + cosTheta);
423  if ( (leftDenominator * rightDenominator) != 0. )
424  {
425  cumul = cumul + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator));
426  }
427  }
428 
429  // Select cosTheta
430  for (G4int i=0; i<iMax; i++)
431  {
432  cosTheta = -1 + i*2./(iMax-1);
433  leftDenominator = (1. + 2.*gamma - cosTheta);
434  rightDenominator = (1. + 2.*delta + cosTheta);
435  if (cumul !=0 && (leftDenominator * rightDenominator) != 0.)
436  value = value + (1./(leftDenominator*leftDenominator) + beta/(rightDenominator*rightDenominator)) / cumul;
437  if (random < value) break;
438  }
439 
440  return cosTheta;
441  }
442 
443  return 0.;
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
447 
449  std::vector<
450  G4double>& vec)
451 {
452  // Sum_{i=0}^{size-1} vector_i k^i
453  //
454  // Phys. Med. Biol. 29 N.4 (1983) 443-447
455 
456  G4double result = 0.;
457  size_t size = vec.size();
458 
459  while (size > 0)
460  {
461  size--;
462 
463  result *= k;
464  result += vec[size];
465  }
466 
467  return result;
468 }
469 
470 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
471 
473  G4double z)
474 {
475 
476  // d sigma_el sigma_Ruth(K)
477  // ------------ (K) ~ -----------------------------
478  // d Omega (1 + 2 n(K) - cos(theta))^2
479  //
480  // We extract cos(theta) distributed as (1 + 2 n(K) - cos(theta))^-2
481  //
482  // 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)
483  //
484  // Phys. Med. Biol. 45 (2000) 3171-3194
485 
486  // ***** Original method
487 
488  if (!fasterCode)
489  {
490 
491  G4double n = ScreeningFactor(k, z);
492 
493  G4double oneOverMax = (4. * n * n);
494 
495  G4double cosTheta = 0.;
496  G4double fCosTheta;
497 
498  do
499  {
500  cosTheta = 2. * G4UniformRand()- 1.;
501  fCosTheta = (1 + 2.*n - cosTheta);
502  if (fCosTheta !=0.) fCosTheta = oneOverMax / (fCosTheta*fCosTheta);
503  }
504  while (fCosTheta < G4UniformRand());
505 
506  return cosTheta;
507  }
508 
509  // ***** Alternative method using cumulative probability
510  if (fasterCode)
511  {
512 
513  G4double cosTheta = -1;
514  G4double cumul = 0;
515  G4double value = 0;
516  G4double n = ScreeningFactor(k, z);
517  G4double fCosTheta;
518 
519  // Number of integration steps in the -1,1 range
520  G4int iMax=200;
521 
522  G4double random = G4UniformRand();
523 
524  // Cumulate differential cross section
525  for (G4int i=0; i<iMax; i++)
526  {
527  cosTheta = -1 + i*2./(iMax-1);
528  fCosTheta = (1 + 2.*n - cosTheta);
529  if (fCosTheta !=0.) cumul = cumul + 1./(fCosTheta*fCosTheta);
530  }
531 
532  // Select cosTheta
533  for (G4int i=0; i<iMax; i++)
534  {
535  cosTheta = -1 + i*2./(iMax-1);
536  fCosTheta = (1 + 2.*n - cosTheta);
537  if (cumul !=0.) value = value + (1./(fCosTheta*fCosTheta)) / cumul;
538  if (random < value) break;
539  }
540  return cosTheta;
541  }
542 
543  return 0.;
544 }
545 
static const double cm
Definition: G4SIunits.hh:118
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
G4double CalculatePolynomial(G4double k, std::vector< G4double > &vec)
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
size_t GetIndex() const
Definition: G4Material.hh:262
G4double RutherfordCrossSection(G4double energy, G4double z)
int G4int
Definition: G4Types.hh:78
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:725
#define G4UniformRand()
Definition: Randomize.hh:97
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:212
static G4DNAMolecularMaterial * Instance()
G4double ScreeningFactor(G4double energy, G4double z)
static const double pi
Definition: G4SIunits.hh:74
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:732
#define DBL_MAX
Definition: templates.hh:83
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
G4double ScreenedRutherfordRandomizeCosTheta(G4double k, G4double z)