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