49namespace Penetration {
53{ -4.06217193e-08, 3.06848412e-06, -9.93217814e-05,
54 1.80172797e-03, -2.01135480e-02, 1.42939448e-01,
55 -6.48348714e-01, 1.85227848e+00, -3.36450378e+00,
56 4.37785068e+00, -4.20557339e+00, 3.81679083e+00,
62{ 7.3144e-05, -2.2474e-03, 3.4555e-02,
63 -4.3574e-01, 2.8954e+00, -1.0381e+00,
69{ 0.2, 0.5, 1, 2, 3, 4, 5, 6, 7,
75{ 17.68*CLHEP::angstrom,
77 28.49*CLHEP::angstrom,
78 45.35*CLHEP::angstrom,
79 70.03*CLHEP::angstrom,
80 98.05*CLHEP::angstrom,
81 120.56*CLHEP::angstrom,
82 132.73*CLHEP::angstrom,
83 142.60*CLHEP::angstrom,
89 137.9*CLHEP::angstrom,
100 for(int8_t i=12; i!=-1 ; --i){
101 r_mean+=
gCoeff[12-i]*std::pow(k_eV,i);
103 r_mean*=CLHEP::nanometer;
114 for(int8_t i=6; i!=-1 ; --i){
115 r_mean+=
gCoeff[6-i]*std::pow(k_eV,i);
117 r_mean*=CLHEP::nanometer;
134 static constexpr double convertRmean3DToSigma1D = 0.62665706865775006;
138 const double sigma1D = r_mean * convertRmean3DToSigma1D;
140 G4RandGauss::shoot(0, sigma1D),
141 G4RandGauss::shoot(0, sigma1D));
170 double r = G4RandGamma::shoot(2,2);
190 return 1e-3*CLHEP::nanometer;
197 description <<
"Terrisol1990 is not tabulated for energies greater than 9eV";
204 size_t lowBin, upBin;
207 lowBin=std::floor(k_eV)+1;
208 upBin=std::min(lowBin+1,
size_t(10));
224 double tanA = (lowS-upS)/(lowE-upE);
225 double sigma3D = lowS + (k_eV-lowE)*tanA;
232 static constexpr double s2r=1.595769121605731;
234 double r_mean=sigma3D*s2r;
242 static constexpr double factor = 2.20496999539;
244 double sigma1D = std::sqrt(std::pow(sigma3D, 2.)*factor);
247 G4RandGauss::shoot(0, sigma1D),
248 G4RandGauss::shoot(0, sigma1D));
258 G4String modelNamePrefix(
"DNAOneStepThermalizationModel_");
260 if(penetrationModel ==
"Terrisol1990")
264 else if(penetrationModel ==
"Meesungnoen2002")
268 else if(penetrationModel ==
"Meesungnoen2002_amorphous")
272 else if(penetrationModel ==
"Kreipl2009")
276 else if(penetrationModel ==
"Ritchie1994")
283 description << penetrationModel +
" is not a valid model name.";
288 "Options are: Terrisol1990, Meesungnoen2002, Ritchie1994.");
301 return Create(
"Ritchie1994");
303 return Create(
"Terrisol1990");
305 return Create(
"Kreipl2009");
307 return Create(
"Meesungnoen2002_amorphous");
310 return Create(
"Meesungnoen2002");
312 G4Exception(
"G4DNASolvationModelFactory::GetMacroDefinedModel",
315 "The solvation parameter stored in G4EmParameters is unknown. Supported types are: fRitchie1994eSolvation, fTerrisol1990eSolvation, fMeesungnoen2002eSolvation.");
@ fMeesungnoen2002eSolvation
@ fMeesungnoensolid2002eSolvation
@ fTerrisol1990eSolvation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
static G4VEmModel * Create(const G4String &penetrationModel)
static G4Electron * Definition()
static G4EmParameters * Instance()
G4DNAModelSubType DNAeSolvationSubType() const
void GetGaussianPenetrationFromRmean3D(G4double r_mean, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double GetRmean(double energy)
static const double gCoeff[7]
static const double gCoeff[13]
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double GetRmean(double energy)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double GetRmean(double energy)
static const double gEnergies_T1990[11]
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static const double gStdDev_T1990[11]
static double Get3DStdDeviation(double energy)