46 fpWaterDensity =
nullptr;
48 lowEnergyLimit = 0 * eV;
49 intermediateEnergyLimit = 200 * eV;
50 highEnergyLimit = 1. * MeV;
66 G4cout <<
"Screened Rutherford Elastic model is constructed "
69 << lowEnergyLimit / eV <<
" eV - "
70 << highEnergyLimit / MeV <<
" MeV"
95 G4cout <<
"Calling G4DNAScreenedRutherfordElasticModel::Initialise()"
102 G4Exception (
"*** WARNING: the G4DNAScreenedRutherfordElasticModel is not "
103 "intented to be used with another particle than the electron",
110 G4Exception(
"*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
111 "not validated below 9 eV",
117 G4Exception(
"*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
118 "not validated above 1 MeV",
126 G4cout <<
"Screened Rutherford elastic model is initialized " <<
G4endl
134 if (isInitialised) {
return; }
141 isInitialised =
true;
201 if (verboseLevel > 3)
203 G4cout <<
"Calling CrossSectionPerVolume() of "
204 "G4DNAScreenedRutherfordElasticModel"
217 G4double n = ScreeningFactor(ekin,z);
218 G4double crossSection = RutherfordCrossSection(ekin, z);
219 sigma = pi * crossSection / (n * (n + 1.));
223 if (verboseLevel > 2)
225 G4cout <<
"__________________________________" <<
G4endl;
226 G4cout <<
"=== G4DNAScreenedRutherfordElasticModel - XS INFO START"
228 G4cout <<
"=== Kinetic energy(eV)=" << ekin/eV
229 <<
" particle : " << particleDefinition->GetParticleName()
231 G4cout <<
"=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
233 G4cout <<
"=== Cross section per water molecule (cm^-1)="
234 << sigma*waterDensity/(1./cm) <<
G4endl;
235 G4cout <<
"=== G4DNAScreenedRutherfordElasticModel - XS INFO END"
240 return sigma*waterDensity;
245G4double G4DNAScreenedRutherfordElasticModel::RutherfordCrossSection(
G4double k,
257 G4double length = (e_squared * (k + electron_mass_c2))
258 / (4 * pi * epsilon0 * k * (k + 2 * electron_mass_c2));
259 G4double cross = z * (z + 1) * length * length;
285 G4double numerator = (alpha_1 + beta_1 *
G4Log(k / eV)) * constK
286 * std::pow(z, 2. / 3.);
288 k /= electron_mass_c2;
293 if (denominator > 0.) value = numerator / denominator;
308 if (verboseLevel > 3)
310 G4cout <<
"Calling SampleSecondaries() of "
311 "G4DNAScreenedRutherfordElasticModel"
319 if (electronEnergy0<intermediateEnergyLimit)
322 if (verboseLevel > 3)
323 {
G4cout <<
"---> Using Brenner & Zaider model" <<
G4endl;}
325 cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
328 if (electronEnergy0>=intermediateEnergyLimit)
331 if (verboseLevel > 3)
332 {
G4cout <<
"---> Using Screened Rutherford model" <<
G4endl;}
335 cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
344 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
346 xDir *= std::cos(phi);
347 yDir *= std::sin(phi);
349 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
358G4double G4DNAScreenedRutherfordElasticModel::
359BrennerZaiderRandomizeCosTheta(
G4double k)
379 gamma = CalculatePolynomial(k, gamma100_200Coeff);
386 gamma =
G4Exp(CalculatePolynomial(k, gamma10_100Coeff));
390 gamma =
G4Exp(CalculatePolynomial(k, gamma035_10Coeff));
400 / (1. / (4. * gamma * gamma) + beta
401 / ((2. + 2. * delta) * (2. + 2. * delta)));
412 leftDenominator = (1. + 2.*gamma - cosTheta);
413 rightDenominator = (1. + 2.*delta + cosTheta);
414 if ( (leftDenominator * rightDenominator) != 0. )
416 fCosTheta = oneOverMax * (1./(leftDenominator*leftDenominator)
417 + beta/(rightDenominator*rightDenominator));
476 G4double S = 2.0 / (a1 * a2) + 2.0 * beta / (b1 * b2);
S = 1.0 /
S;
479 G4double A =
S * (b1 - beta * a2) + cp * a2 * b1;
480 G4double B =
S * (b1 * b2 + beta * a1 * a2) - cp * a2 * b1 * c1;
481 G4double C =
S * (b * b1 + beta * a * a2) - cp * a2 * b1 * c2;
484 return (-1.0 * B + std::sqrt(B * B - 4.0 *
A * C)) / (2.0 *
A);
530G4double G4DNAScreenedRutherfordElasticModel::
532 std::vector<G4double>& vec)
539 size_t size = vec.size();
554G4double G4DNAScreenedRutherfordElasticModel::
555ScreenedRutherfordRandomizeCosTheta(
G4double k,
584 fCosTheta = (1 + 2.*
n - cosTheta);
585 if (fCosTheta !=0.) fCosTheta = oneOverMax / (fCosTheta*fCosTheta);
616 G4double numerator = cp * (1.0 + 2.0 *
n) - n;
618 return numerator / denominator;
G4double C(G4double temp)
G4double S(G4double temp)
G4double B(G4double temperature)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
G4GLOB_DLL std::ostream G4cout
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
static G4DNAMolecularMaterial * Instance()
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ParticleChangeForGamma * fParticleChangeForGamma
~G4DNAScreenedRutherfordElasticModel() override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4DNAScreenedRutherfordElasticModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAScreenedRutherfordElasticModel")
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
std::size_t GetIndex() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)