source: trunk/source/materials/src/G4IonisParamElm.cc@ 1215

Last change on this file since 1215 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 8.7 KB
Line 
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,v 1.15 2008/06/03 14:30:10 vnivanch Exp $
28// GEANT4 tag $Name: materials-V09-02-18 $
29//
30//
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
32//
33// 06-09-04, Update calculated values after any change of ionisation
34// potential change. V.Ivanchenko
35// 17-10-02, Fix excitation energy interpolation. V.Ivanchenko
36// 08-03-01, correct handling of fShellCorrectionVector. M.Maire
37// 22-11-00, tabulation of ionisation potential from
38// the ICRU Report N#37. V.Ivanchenko
39// 09-07-98, data moved from G4Element. M.Maire
40//
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
42
43#include "G4IonisParamElm.hh"
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
46
47G4IonisParamElm::G4IonisParamElm(G4double Z)
48{
49 if (Z < 1.)
50 G4Exception (" ERROR! It is not allowed to create an Element with Z < 1" );
51
52 // some basic functions of the atomic number
53 fZ = Z;
54 fZ3 = std::pow(fZ, 1./3.);
55 fZZ3 = std::pow(fZ*(fZ+1.), 1./3.);
56 flogZ3 = std::log(fZ)/3.;
57
58 // Parameters for energy loss by ionisation
59
60 // Mean excitation energy
61 // from "Stopping Powers for Electrons and Positrons"
62 // ICRU Report N#37, 1984 (energy in eV)
63 static double exc[100] = {
64 21.8, 20.9, 13.3, 15.9, 15.2, 13.0, 11.7, 11.9, 11.5, 13.7,
65 13.6, 13.0, 12.8, 12.4, 11.5, 11.3, 10.2, 10.4, 10.0, 9.6,
66 10.3, 10.6, 10.7, 10.7, 10.9, 11.0, 11.0, 11.1, 11.1, 11.0,
67 10.8, 10.4, 10.5, 10.2, 9.8, 9.8, 9.8, 9.6, 9.7, 9.8,
68 10.2, 10.1, 10.0, 10.0, 10.0, 10.2, 10.0, 9.8, 10.0, 9.8,
69 9.5, 9.3, 9.3, 8.9, 8.9, 8.8, 8.8, 8.8, 9.1, 9.1,
70 9.2, 9.3, 9.2, 9.2, 9.4, 9.5, 9.7, 9.7, 9.8, 9.8,
71 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 10.1, 10.0, 10.0,
72 10.0, 10.0, 9.9, 9.9, 9.7, 9.2, 9.5, 9.4, 9.4, 9.4,
73 9.6, 9.7, 9.7, 9.8, 9.8, 9.8, 9.8, 9.9, 9.9, 9.9 };
74
75 G4int iz = (G4int)Z - 1 ;
76 if(0 > iz) iz = 0;
77 else if(99 < iz) iz = 99 ;
78 fMeanExcitationEnergy = fZ * exc[iz] * eV ;
79
80 // compute parameters for ion transport
81 // The aproximation from:
82 // J.F.Ziegler, J.P. Biersack, U. Littmark
83 // The Stopping and Range of Ions in Matter,
84 // Vol.1, Pergamon Press, 1985
85 // Fast ions or hadrons
86
87 if(91 < iz) iz = 91;
88
89 static G4double vFermi[92] = {
90 1.0309, 0.15976, 0.59782, 1.0781, 1.0486, 1.0, 1.058, 0.93942, 0.74562, 0.3424,
91 0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712,
92 0.81707, 0.9943, 1.1423, 1.2381, 1.1222, 0.92705, 1.0047, 1.2, 1.0661, 0.97411,
93 0.84912, 0.95, 1.0903, 1.0429, 0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207,
94 1.029, 1.2542, 1.122, 1.1241, 1.0882, 1.2709, 1.2542, 0.90094, 0.74093, 0.86054,
95 0.93155, 1.0047, 0.55379, 0.43289, 0.32636, 0.5131, 0.695, 0.72591, 0.71202, 0.67413,
96 0.71418, 0.71453, 0.5911, 0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884,
97 0.71801, 0.83048, 1.1222, 1.2381, 1.045, 1.0733, 1.0953, 1.2381, 1.2879, 0.78654,
98 0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065,
99 1.9578, 1.0257} ;
100
101 static G4double lFactor[92] = {
102 1.0, 1.0, 1.1, 1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9,
103 0.82, 0.81, 0.83, 0.88, 1.0, 0.95, 0.97, 0.99, 0.98, 0.97,
104 0.98, 0.97, 0.96, 0.93, 0.91, 0.9, 0.88, 0.9, 0.9, 0.9,
105 0.9, 0.85, 0.9, 0.9, 0.91, 0.92, 0.9, 0.9, 0.9, 0.9,
106 0.9, 0.88, 0.9, 0.88, 0.88, 0.9, 0.9, 0.88, 0.9, 0.9,
107 0.9, 0.9, 0.96, 1.2, 0.9, 0.88, 0.88, 0.85, 0.9, 0.9,
108 0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1, 1.08, 1.08,
109 1.08, 1.08, 1.09, 1.09, 1.1, 1.11, 1.12, 1.13, 1.14, 1.15,
110 1.17, 1.2, 1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16,
111 1.16, 1.16} ;
112
113 fVFermi = vFermi[iz];
114 fLFactor = lFactor[iz];
115
116 // obsolete parameters for ionisation
117 fTau0 = 0.1*fZ3*MeV/proton_mass_c2;
118 fTaul = 2.*MeV/proton_mass_c2;
119
120 // compute the Bethe-Bloch formula for energy = fTaul*particle mass
121 G4double rate = fMeanExcitationEnergy/electron_mass_c2 ;
122 G4double w = fTaul*(fTaul+2.) ;
123 fBetheBlochLow = (fTaul+1.)*(fTaul+1.)*std::log(2.*w/rate)/w - 1. ;
124 fBetheBlochLow = 2.*fZ*twopi_mc2_rcl2*fBetheBlochLow ;
125
126 fClow = std::sqrt(fTaul)*fBetheBlochLow;
127 fAlow = 6.458040 * fClow/fTau0;
128 G4double Taum = 0.035*fZ3*MeV/proton_mass_c2;
129 fBlow =-3.229020*fClow/(fTau0*std::sqrt(Taum));
130
131 // Shell correction parameterization
132 fShellCorrectionVector = new G4double[3]; //[3]
133 rate = 0.001*fMeanExcitationEnergy/eV;
134 G4double rate2 = rate*rate;
135 /*
136 fShellCorrectionVector[0] = ( 1.10289e5 + 5.14781e8*rate)*rate2 ;
137 fShellCorrectionVector[1] = ( 7.93805e3 - 2.22565e7*rate)*rate2 ;
138 fShellCorrectionVector[2] = (-9.92256e1 + 2.10823e5*rate)*rate2 ;
139 */
140 fShellCorrectionVector[0] = ( 0.422377 + 3.858019*rate)*rate2 ;
141 fShellCorrectionVector[1] = ( 0.0304043 - 0.1667989*rate)*rate2 ;
142 fShellCorrectionVector[2] = (-0.00038106 + 0.00157955*rate)*rate2 ;
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
146
147// Fake default constructor - sets only member data and allocates memory
148// for usage restricted to object persistency
149
150G4IonisParamElm::G4IonisParamElm(__void__&)
151 : fShellCorrectionVector(0)
152{
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
156
157G4IonisParamElm::~G4IonisParamElm()
158{
159 if (fShellCorrectionVector) delete [] fShellCorrectionVector;
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
163
164G4IonisParamElm::G4IonisParamElm(G4IonisParamElm& right)
165{
166 *this = right;
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
170
171const G4IonisParamElm& G4IonisParamElm::operator=(const G4IonisParamElm& right)
172{
173 if (this != &right)
174 {
175 fZ = right.fZ;
176 fZ3 = right.fZ3;
177 fZZ3 = right.fZZ3;
178 flogZ3 = right.flogZ3;
179 fTau0 = right.fTau0;
180 fTaul = right.fTaul;
181 fBetheBlochLow = right.fBetheBlochLow;
182 fAlow = right.fAlow;
183 fBlow = right.fBlow;
184 fClow = right.fClow;
185 fMeanExcitationEnergy = right.fMeanExcitationEnergy;
186 if (fShellCorrectionVector) delete [] fShellCorrectionVector;
187 fShellCorrectionVector = new G4double[3];
188 fShellCorrectionVector[0] = right.fShellCorrectionVector[0];
189 fShellCorrectionVector[1] = right.fShellCorrectionVector[1];
190 fShellCorrectionVector[2] = right.fShellCorrectionVector[2];
191 }
192 return *this;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
196
197G4int G4IonisParamElm::operator==(const G4IonisParamElm& right) const
198{
199 return (this == (G4IonisParamElm *) &right);
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
203
204G4int G4IonisParamElm::operator!=(const G4IonisParamElm& right) const
205{
206 return (this != (G4IonisParamElm *) &right);
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.