source: trunk/source/processes/electromagnetic/utils/src/G4ionEffectiveCharge.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 7.5 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// $Id: G4ionEffectiveCharge.cc,v 1.24 2008/12/18 13:01:46 gunter Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4ionEffectiveCharge
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 07.05.2002
39//
40// Modifications:
41// 12.09.2004 Set low energy limit to 1 keV (V.Ivanchenko)
42// 25.01.2005 Add protection - min Charge 0.1 eplus (V.Ivanchenko)
43// 28.04.2006 Set upper energy limit to 50 MeV (V.Ivanchenko)
44// 23.05.2006 Set upper energy limit to Z*10 MeV (V.Ivanchenko)
45// 15.08.2006 Add protection for not defined material (V.Ivanchenko)
46// 27-09-2007 Use Fermi energy from material, optimazed formulas (V.Ivanchenko)
47//
48
49// -------------------------------------------------------------------
50//
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
54#include "G4ionEffectiveCharge.hh"
55#include "G4UnitsTable.hh"
56#include "G4Material.hh"
57#include "G4NistManager.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60
61G4ionEffectiveCharge::G4ionEffectiveCharge()
62{
63 chargeCorrection = 1.0;
64 energyHighLimit = 20.0*MeV;
65 energyLowLimit = 1.0*keV;
66 energyBohr = 25.*keV;
67 massFactor = amu_c2/(proton_mass_c2*keV);
68 minCharge = 1.0;
69 lastPart = 0;
70 lastMat = 0;
71 lastKinEnergy = 0.0;
72 effCharge = eplus;
73 nist = G4NistManager::Instance();
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78G4ionEffectiveCharge::~G4ionEffectiveCharge()
79{}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
83G4double G4ionEffectiveCharge::EffectiveCharge(const G4ParticleDefinition* p,
84 const G4Material* material,
85 G4double kineticEnergy)
86{
87 if(p == lastPart && material == lastMat && kineticEnergy == lastKinEnergy)
88 return effCharge;
89
90 lastPart = p;
91 lastMat = material;
92 lastKinEnergy = kineticEnergy;
93
94 G4double mass = p->GetPDGMass();
95 G4double charge = p->GetPDGCharge();
96 G4double Zi = charge/eplus;
97
98 chargeCorrection = 1.0;
99 effCharge = charge;
100
101 // The aproximation of ion effective charge from:
102 // J.F.Ziegler, J.P. Biersack, U. Littmark
103 // The Stopping and Range of Ions in Matter,
104 // Vol.1, Pergamon Press, 1985
105 // Fast ions or hadrons
106 G4double reducedEnergy = kineticEnergy * proton_mass_c2/mass ;
107 if( reducedEnergy > Zi*energyHighLimit || Zi < 1.5 || !material) return charge;
108
109 G4double z = material->GetIonisation()->GetZeffective();
110 // reducedEnergy = std::max(reducedEnergy,energyLowLimit);
111
112 // Helium ion case
113 if( Zi < 2.5 ) {
114
115 static G4double c[6] = {0.2865, 0.1266, -0.001429,
116 0.02402,-0.01135, 0.001475} ;
117
118 G4double Q = std::max(0.0,std::log(reducedEnergy*massFactor));
119 G4double x = c[0];
120 G4double y = 1.0;
121 for (G4int i=1; i<6; i++) {
122 y *= Q;
123 x += y * c[i] ;
124 }
125 G4double ex;
126 if(x < 0.2) ex = x * (1 - 0.5*x);
127 else ex = 1. - std::exp(-x);
128
129 G4double tq = 7.6 - Q;
130 G4double tq2= tq*tq;
131 G4double tt = ( 0.007 + 0.00005 * z );
132 if(tq2 < 0.2) tt *= (1.0 - tq2 + 0.5*tq2*tq2);
133 else tt *= std::exp(-tq2);
134
135 effCharge = charge*(1.0 + tt) * std::sqrt(ex);
136
137 // Heavy ion case
138 } else {
139
140 G4double y;
141 // = nist->GetZ13(z);
142 //G4double z23 = y*y;
143 G4double zi13 = nist->GetZ13(Zi);
144 G4double zi23 = zi13*zi13;
145 // G4double e = std::max(reducedEnergy,energyBohr/z23);
146 //G4double e = reducedEnergy;
147
148 // v1 is ion velocity in vF unit
149 G4double eF = material->GetIonisation()->GetFermiEnergy();
150 G4double v1sq = reducedEnergy/eF;
151 G4double vFsq = eF/energyBohr;
152 G4double vF = std::sqrt(eF/energyBohr);
153
154 // Faster than Fermi velocity
155 if ( v1sq > 1.0 ) {
156 y = vF * std::sqrt(v1sq) * ( 1.0 + 0.2/v1sq ) / zi23 ;
157
158 // Slower than Fermi velocity
159 } else {
160 y = 0.692308 * vF * (1.0 + 0.666666*v1sq + v1sq*v1sq/15.0) / zi23 ;
161 }
162
163 G4double q;
164 G4double y3 = std::pow(y, 0.3) ;
165 // G4cout<<"y= "<<y<<" y3= "<<y3<<" v1= "<<v1<<" vF= "<<vF<<G4endl;
166 q = 1.0 - std::exp( 0.803*y3 - 1.3167*y3*y3 - 0.38157*y - 0.008983*y*y ) ;
167
168 //y *= 0.77;
169 //y *= (0.75 + 0.52/Zi);
170
171 //if( y < 0.2 ) q = y*(1.0 - 0.5*y);
172 //else q = 1.0 - std::exp(-y);
173
174 G4double qmin = minCharge/Zi;
175 if(q < qmin) q = qmin;
176
177 effCharge = q*charge;
178
179 /*
180 G4double x1 = 1.0*effCharge*(1.0 - 0.132*std::log(y))/(y*std::sqrt(z));
181 G4double x2 = 0.1*effCharge*effCharge*energyBohr/reducedEnergy;
182
183 chargeCorrection = 1.0 + x1 - x2;
184
185 G4cout << "x1= "<<x1<<" x2= "<< x2<<" corr= "<<chargeCorrection<<G4endl;
186 */
187
188 G4double tq = 7.6 - std::log(reducedEnergy/keV);
189 G4double tq2= tq*tq;
190 G4double sq = ( 0.18 + 0.0015 * z ) / (Zi*Zi);
191 if(tq2 < 0.2) sq *= (1.0 - tq2 + 0.5*tq2*tq2);
192 else sq *= std::exp(-tq2);
193 sq += 1.0;
194 // G4cout << "sq= " << sq << G4endl;
195
196 // Screen length according to
197 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
198 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
199
200 G4double lambda = 10.0 * vF / (zi13 * (6.0 + q));
201 if(q < 0.2) lambda *= (1.0 - 0.66666667*q - q*q/9.0);
202 else lambda *= std::pow(1.0-q, 0.666666);
203
204 G4double lambda2 = lambda*lambda;
205
206 G4double xx = (0.5/q - 0.5)/vFsq;
207 if(lambda2 < 0.2) xx *= lambda2*(1.0 - 0.5*lambda2);
208 else xx *= std::log(1.0 + lambda2);
209
210 chargeCorrection = sq * (1.0 + xx);
211
212 }
213 // G4cout << "G4ionEffectiveCharge: charge= " << charge << " q= " << q
214 // << " chargeCor= " << chargeCorrection
215 // << " e(MeV)= " << kineticEnergy/MeV << G4endl;
216 return effCharge;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220
221
Note: See TracBrowser for help on using the repository browser.