Geant4  10.00.p02
G4IonisParamElm.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 //
27 // $Id: G4IonisParamElm.cc 68721 2013-04-05 09:20:30Z gcosmo $
28 //
29 //
30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
31 //
32 // 09-07-98, data moved from G4Element. M.Maire
33 // 22-11-00, tabulation of ionisation potential from
34 // the ICRU Report N#37. V.Ivanchenko
35 // 08-03-01, correct handling of fShellCorrectionVector. M.Maire
36 // 17-10-02, Fix excitation energy interpolation. V.Ivanchenko
37 // 06-09-04, Update calculated values after any change of ionisation
38 // potential change. V.Ivanchenko
39 // 29-04-10, Using G4Pow and mean ionisation energy from NIST V.Ivanchenko
40 // 27.10.11: new scheme for G4Exception (mma)
41 //
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
43 
44 #include "G4IonisParamElm.hh"
45 #include "G4NistManager.hh"
46 #include "G4Pow.hh"
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
51 
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 
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 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
137 
138 // Fake default constructor - sets only member data and allocates memory
139 // for usage restricted to object persistency
140 
142  : fShellCorrectionVector(0)
143 {
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
149 
151 {
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
156 
158 {
160  *this = right;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
164 
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;
176  fAlow = right.fAlow;
177  fBlow = right.fBlow;
178  fClow = right.fClow;
180  fVFermi = right.fVFermi;
181  fLFactor = right.fLFactor;
187  }
188  return *this;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
192 
194 {
195  return (this == (G4IonisParamElm *) &right);
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
199 
201 {
202  return (this != (G4IonisParamElm *) &right);
203 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
static const double MeV
Definition: G4SIunits.hh:193
virtual ~G4IonisParamElm()
Definition: G4Pow.hh:56
int G4int
Definition: G4Types.hh:78
G4int operator==(const G4IonisParamElm &) const
static G4NistManager * Instance()
G4double logZ(G4int Z) const
Definition: G4Pow.hh:165
G4double Z13(G4int Z) const
Definition: G4Pow.hh:129
G4double fMeanExcitationEnergy
G4double iz
Definition: TRTMaterials.hh:39
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:194
G4IonisParamElm(G4double Z)
G4int operator!=(const G4IonisParamElm &) const
G4IonisParamElm & operator=(const G4IonisParamElm &)
G4double * fShellCorrectionVector
double G4double
Definition: G4Types.hh:76
G4double GetMeanIonisationEnergy(G4int Z) const