Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4IonisParamElm Class Reference

#include <G4IonisParamElm.hh>

Public Member Functions

 G4IonisParamElm (G4double Z)
 
virtual ~G4IonisParamElm ()
 
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
 
 G4IonisParamElm (G4IonisParamElm &)
 
G4IonisParamElmoperator= (const G4IonisParamElm &)
 
G4int operator== (const G4IonisParamElm &) const
 
G4int operator!= (const G4IonisParamElm &) const
 
 G4IonisParamElm (__void__ &)
 

Detailed Description

Definition at line 51 of file G4IonisParamElm.hh.

Constructor & Destructor Documentation

G4IonisParamElm::G4IonisParamElm ( G4double  Z)

Definition at line 52 of file G4IonisParamElm.cc.

References python.hepunit::electron_mass_c2, python.hepunit::eV, FatalException, G4Exception(), G4Pow::GetInstance(), G4NistManager::GetMeanIonisationEnergy(), G4NistManager::Instance(), iz, G4Pow::logZ(), python.hepunit::MeV, python.hepunit::proton_mass_c2, python.hepunit::twopi_mc2_rcl2, and G4Pow::Z13().

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

Definition at line 150 of file G4IonisParamElm.cc.

151 {
152  if (fShellCorrectionVector) { delete [] fShellCorrectionVector; }
153 }
G4IonisParamElm::G4IonisParamElm ( G4IonisParamElm right)

Definition at line 157 of file G4IonisParamElm.cc.

References right.

158 {
159  fShellCorrectionVector = 0;
160  *this = right;
161 }
G4IonisParamElm::G4IonisParamElm ( __void__ &  )

Definition at line 141 of file G4IonisParamElm.cc.

142  : fShellCorrectionVector(0)
143 {
144  fZ=fZ3=fZZ3=flogZ3=fTau0=fTaul=fBetheBlochLow=fAlow=fBlow=fClow
145  =fMeanExcitationEnergy=fVFermi=fLFactor=0.0;
146 }

Member Function Documentation

G4double G4IonisParamElm::GetAlow ( ) const
inline

Definition at line 70 of file G4IonisParamElm.hh.

70 {return fAlow;}
G4double G4IonisParamElm::GetBlow ( ) const
inline

Definition at line 71 of file G4IonisParamElm.hh.

71 {return fBlow;}
G4double G4IonisParamElm::GetClow ( ) const
inline

Definition at line 72 of file G4IonisParamElm.hh.

72 {return fClow;}
G4double G4IonisParamElm::GetFermiVelocity ( ) const
inline

Definition at line 78 of file G4IonisParamElm.hh.

78 {return fVFermi;};
G4double G4IonisParamElm::GetLFactor ( ) const
inline

Definition at line 79 of file G4IonisParamElm.hh.

79 {return fLFactor;};
G4double G4IonisParamElm::GetlogZ3 ( ) const
inline
G4double G4IonisParamElm::GetMeanExcitationEnergy ( ) const
inline

Definition at line 75 of file G4IonisParamElm.hh.

75 {return fMeanExcitationEnergy;}
G4double* G4IonisParamElm::GetShellCorrectionVector ( ) const
inline

Definition at line 81 of file G4IonisParamElm.hh.

81 {return fShellCorrectionVector;}
G4double G4IonisParamElm::GetTau0 ( ) const
inline

Definition at line 65 of file G4IonisParamElm.hh.

65 {return fTau0;};
G4double G4IonisParamElm::GetTaul ( ) const
inline

Definition at line 68 of file G4IonisParamElm.hh.

68 {return fTaul;} // 2*MeV/proton mass
G4double G4IonisParamElm::GetZ ( ) const
inline

Definition at line 60 of file G4IonisParamElm.hh.

60 {return fZ;} // atomic number Z
G4double G4IonisParamElm::GetZ3 ( ) const
inline
G4double G4IonisParamElm::GetZZ3 ( ) const
inline

Definition at line 62 of file G4IonisParamElm.hh.

62 {return fZZ3;} // std::pow (Z(Z+1),1/3)
G4int G4IonisParamElm::operator!= ( const G4IonisParamElm right) const

Definition at line 200 of file G4IonisParamElm.cc.

201 {
202  return (this != (G4IonisParamElm *) &right);
203 }
G4IonisParamElm & G4IonisParamElm::operator= ( const G4IonisParamElm right)

Definition at line 165 of file G4IonisParamElm.cc.

166 {
167  if (this != &right)
168  {
169  fZ = right.fZ;
170  fZ3 = right.fZ3;
171  fZZ3 = right.fZZ3;
172  flogZ3 = right.flogZ3;
173  fTau0 = right.fTau0;
174  fTaul = right.fTaul;
175  fBetheBlochLow = right.fBetheBlochLow;
176  fAlow = right.fAlow;
177  fBlow = right.fBlow;
178  fClow = right.fClow;
179  fMeanExcitationEnergy = right.fMeanExcitationEnergy;
180  fVFermi = right.fVFermi;
181  fLFactor = right.fLFactor;
182  if (fShellCorrectionVector) { delete [] fShellCorrectionVector; }
183  fShellCorrectionVector = new G4double[3];
184  fShellCorrectionVector[0] = right.fShellCorrectionVector[0];
185  fShellCorrectionVector[1] = right.fShellCorrectionVector[1];
186  fShellCorrectionVector[2] = right.fShellCorrectionVector[2];
187  }
188  return *this;
189 }
double G4double
Definition: G4Types.hh:76
G4int G4IonisParamElm::operator== ( const G4IonisParamElm right) const

Definition at line 193 of file G4IonisParamElm.cc.

194 {
195  return (this == (G4IonisParamElm *) &right);
196 }

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