Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IonisParamElm Class Reference

#include <G4IonisParamElm.hh>

Public Member Functions

 G4IonisParamElm (G4double Z)
 
 ~G4IonisParamElm ()
 
G4IonisParamElmoperator= (const G4IonisParamElm &)=delete
 
 G4IonisParamElm (G4IonisParamElm &)=delete
 
G4double GetZ () const
 
G4double GetZ3 () const
 
G4double GetZZ3 () const
 
G4double GetlogZ3 () const
 
G4double GetTau0 () const
 
G4double GetTaul () const
 
G4double GetAlow () const
 
G4double GetBlow () const
 
G4double GetClow () const
 
G4double GetMeanExcitationEnergy () const
 
G4double GetFermiVelocity () const
 
G4double GetLFactor () const
 
G4doubleGetShellCorrectionVector () const
 
G4bool operator== (const G4IonisParamElm &) const =delete
 
G4bool operator!= (const G4IonisParamElm &) const =delete
 

Detailed Description

Definition at line 40 of file G4IonisParamElm.hh.

Constructor & Destructor Documentation

◆ G4IonisParamElm() [1/2]

G4IonisParamElm::G4IonisParamElm ( G4double Z)

Definition at line 45 of file G4IonisParamElm.cc.

46{
47 G4int Z = G4lrint(AtomNumber);
48 if (Z < 1) {
49 G4Exception("G4IonisParamElm::G4IonisParamElm()", "mat501", FatalException,
50 "It is not allowed to create an Element with Z<1");
51 }
52 G4Pow* g4pow = G4Pow::GetInstance();
53
54 // some basic functions of the atomic number
55 fZ = Z;
56 fZ3 = g4pow->Z13(Z);
57 fZZ3 = fZ3 * g4pow->Z13(Z + 1);
58 flogZ3 = g4pow->logZ(Z) / 3.;
59
60 fMeanExcitationEnergy = G4NistManager::Instance()->GetMeanIonisationEnergy(Z);
61
62 // compute parameters for ion transport
63 // The aproximation from:
64 // J.F.Ziegler, J.P. Biersack, U. Littmark
65 // The Stopping and Range of Ions in Matter,
66 // Vol.1, Pergamon Press, 1985
67 // Fast ions or hadrons
68
69 G4int iz = Z - 1;
70 if (91 < iz) {
71 iz = 91;
72 }
73
74 // clang-format off
75 static const G4double vFermi[92] = {
76 1.0309, 0.15976, 0.59782, 1.0781, 1.0486, 1.0, 1.058, 0.93942, 0.74562, 0.3424,
77 0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712,
78 0.81707, 0.9943, 1.1423, 1.2381, 1.1222, 0.92705, 1.0047, 1.2, 1.0661, 0.97411,
79 0.84912, 0.95, 1.0903, 1.0429, 0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207,
80 1.029, 1.2542, 1.122, 1.1241, 1.0882, 1.2709, 1.2542, 0.90094, 0.74093, 0.86054,
81 0.93155, 1.0047, 0.55379, 0.43289, 0.32636, 0.5131, 0.695, 0.72591, 0.71202, 0.67413,
82 0.71418, 0.71453, 0.5911, 0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884,
83 0.71801, 0.83048, 1.1222, 1.2381, 1.045, 1.0733, 1.0953, 1.2381, 1.2879, 0.78654,
84 0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065,
85 1.9578, 1.0257} ;
86
87 static const G4double lFactor[92] = {
88 1.0, 1.0, 1.1, 1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9,
89 0.82, 0.81, 0.83, 0.88, 1.0, 0.95, 0.97, 0.99, 0.98, 0.97,
90 0.98, 0.97, 0.96, 0.93, 0.91, 0.9, 0.88, 0.9, 0.9, 0.9,
91 0.9, 0.85, 0.9, 0.9, 0.91, 0.92, 0.9, 0.9, 0.9, 0.9,
92 0.9, 0.88, 0.9, 0.88, 0.88, 0.9, 0.9, 0.88, 0.9, 0.9,
93 0.9, 0.9, 0.96, 1.2, 0.9, 0.88, 0.88, 0.85, 0.9, 0.9,
94 0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1, 1.08, 1.08,
95 1.08, 1.08, 1.09, 1.09, 1.1, 1.11, 1.12, 1.13, 1.14, 1.15,
96 1.17, 1.2, 1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16,
97 1.16, 1.16} ;
98 // clang-format on
99
100 fVFermi = vFermi[iz];
101 fLFactor = lFactor[iz];
102
103 // obsolete parameters for ionisation
104 fTau0 = 0.1 * fZ3 * MeV / proton_mass_c2;
105 fTaul = 2. * MeV / proton_mass_c2;
106
107 // compute the Bethe-Bloch formula for energy = fTaul*particle mass
108 G4double rate = fMeanExcitationEnergy / electron_mass_c2;
109 G4double w = fTaul * (fTaul + 2.);
110 fBetheBlochLow = (fTaul + 1.) * (fTaul + 1.) * std::log(2. * w / rate) / w - 1.;
111 fBetheBlochLow = 2. * fZ * twopi_mc2_rcl2 * fBetheBlochLow;
112
113 fClow = std::sqrt(fTaul) * fBetheBlochLow;
114 fAlow = 6.458040 * fClow / fTau0;
115 G4double Taum = 0.035 * fZ3 * MeV / proton_mass_c2;
116 fBlow = -3.229020 * fClow / (fTau0 * std::sqrt(Taum));
117
118 // Shell correction parameterization
119 fShellCorrectionVector = new G4double[3]; //[3]
120 rate = 0.001 * fMeanExcitationEnergy / eV;
121 G4double rate2 = rate * rate;
122 fShellCorrectionVector[0] = (0.422377 + 3.858019 * rate) * rate2;
123 fShellCorrectionVector[1] = (0.0304043 - 0.1667989 * rate) * rate2;
124 fShellCorrectionVector[2] = (-0.00038106 + 0.00157955 * rate) * rate2;
125}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetMeanIonisationEnergy(G4int Z) const
static G4NistManager * Instance()
Definition G4Pow.hh:49
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double logZ(G4int Z) const
Definition G4Pow.hh:137
G4double Z13(G4int Z) const
Definition G4Pow.hh:123
int G4lrint(double ad)
Definition templates.hh:134

◆ ~G4IonisParamElm()

G4IonisParamElm::~G4IonisParamElm ( )

Definition at line 129 of file G4IonisParamElm.cc.

129{ delete[] fShellCorrectionVector; }

◆ G4IonisParamElm() [2/2]

G4IonisParamElm::G4IonisParamElm ( G4IonisParamElm & )
delete

Member Function Documentation

◆ GetAlow()

G4double G4IonisParamElm::GetAlow ( ) const
inline

Definition at line 69 of file G4IonisParamElm.hh.

69{ return fAlow; }

◆ GetBlow()

G4double G4IonisParamElm::GetBlow ( ) const
inline

Definition at line 70 of file G4IonisParamElm.hh.

70{ return fBlow; }

◆ GetClow()

G4double G4IonisParamElm::GetClow ( ) const
inline

Definition at line 71 of file G4IonisParamElm.hh.

71{ return fClow; }

◆ GetFermiVelocity()

G4double G4IonisParamElm::GetFermiVelocity ( ) const
inline

Definition at line 76 of file G4IonisParamElm.hh.

76{ return fVFermi; };

◆ GetLFactor()

G4double G4IonisParamElm::GetLFactor ( ) const
inline

Definition at line 77 of file G4IonisParamElm.hh.

77{ return fLFactor; };

◆ GetlogZ3()

◆ GetMeanExcitationEnergy()

G4double G4IonisParamElm::GetMeanExcitationEnergy ( ) const
inline

Definition at line 74 of file G4IonisParamElm.hh.

74{ return fMeanExcitationEnergy; }

Referenced by G4ChannelingFastSimCrystalData::SetMaterialProperties().

◆ GetShellCorrectionVector()

G4double * G4IonisParamElm::GetShellCorrectionVector ( ) const
inline

Definition at line 80 of file G4IonisParamElm.hh.

80{ return fShellCorrectionVector; }

◆ GetTau0()

G4double G4IonisParamElm::GetTau0 ( ) const
inline

Definition at line 63 of file G4IonisParamElm.hh.

63{ return fTau0; };

◆ GetTaul()

G4double G4IonisParamElm::GetTaul ( ) const
inline

Definition at line 66 of file G4IonisParamElm.hh.

66{ return fTaul; }

◆ GetZ()

G4double G4IonisParamElm::GetZ ( ) const
inline

Definition at line 51 of file G4IonisParamElm.hh.

51{ return fZ; }

◆ GetZ3()

G4double G4IonisParamElm::GetZ3 ( ) const
inline

Definition at line 54 of file G4IonisParamElm.hh.

54{ return fZ3; }

Referenced by G4BetheHeitler5DModel::SampleSecondaries(), and G4BetheHeitlerModel::SampleSecondaries().

◆ GetZZ3()

G4double G4IonisParamElm::GetZZ3 ( ) const
inline

Definition at line 57 of file G4IonisParamElm.hh.

57{ return fZZ3; }

◆ operator!=()

G4bool G4IonisParamElm::operator!= ( const G4IonisParamElm & ) const
delete

◆ operator=()

G4IonisParamElm & G4IonisParamElm::operator= ( const G4IonisParamElm & )
delete

◆ operator==()

G4bool G4IonisParamElm::operator== ( const G4IonisParamElm & ) const
delete

The documentation for this class was generated from the following files: