Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ionEffectiveCharge.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 // $Id: G4ionEffectiveCharge.cc 74376 2013-10-04 08:25:47Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4ionEffectiveCharge
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 07.05.2002
38 //
39 // Modifications:
40 // 12.09.2004 Set low energy limit to 1 keV (V.Ivanchenko)
41 // 25.01.2005 Add protection - min Charge 0.1 eplus (V.Ivanchenko)
42 // 28.04.2006 Set upper energy limit to 50 MeV (V.Ivanchenko)
43 // 23.05.2006 Set upper energy limit to Z*10 MeV (V.Ivanchenko)
44 // 15.08.2006 Add protection for not defined material (V.Ivanchenko)
45 // 27-09-2007 Use Fermi energy from material, optimazed formulas (V.Ivanchenko)
46 //
47 
48 // -------------------------------------------------------------------
49 //
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 
53 #include "G4ionEffectiveCharge.hh"
54 #include "G4PhysicalConstants.hh"
55 #include "G4SystemOfUnits.hh"
56 #include "G4UnitsTable.hh"
57 #include "G4Material.hh"
58 #include "G4NistManager.hh"
59 #include "G4Log.hh"
60 #include "G4Exp.hh"
61 #include "G4Pow.hh"
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64 
66 {
67  chargeCorrection = 1.0;
68  energyHighLimit = 20.0*MeV;
69  energyLowLimit = 1.0*keV;
70  energyBohr = 25.*keV;
71  massFactor = amu_c2/(proton_mass_c2*keV);
72  minCharge = 1.0;
73  lastPart = 0;
74  lastMat = 0;
75  lastKinEnergy = 0.0;
76  effCharge = eplus;
77  g4pow = G4Pow::GetInstance();
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83 {}
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88  const G4Material* material,
89  G4double kineticEnergy)
90 {
91  if(p == lastPart && material == lastMat && kineticEnergy == lastKinEnergy)
92  return effCharge;
93 
94  lastPart = p;
95  lastMat = material;
96  lastKinEnergy = kineticEnergy;
97 
98  G4double mass = p->GetPDGMass();
99  G4double charge = p->GetPDGCharge();
100  G4double Zi = charge/eplus;
101 
102  chargeCorrection = 1.0;
103  effCharge = charge;
104 
105  // The aproximation of ion effective charge from:
106  // J.F.Ziegler, J.P. Biersack, U. Littmark
107  // The Stopping and Range of Ions in Matter,
108  // Vol.1, Pergamon Press, 1985
109  // Fast ions or hadrons
110  G4double reducedEnergy = kineticEnergy * proton_mass_c2/mass ;
111 
112  //G4cout << "e= " << reducedEnergy << " Zi= " << Zi << " " << material->GetName() << G4endl;
113 
114  if( reducedEnergy > Zi*energyHighLimit || Zi < 1.5 || !material) return charge;
115 
116  G4double z = material->GetIonisation()->GetZeffective();
117  reducedEnergy = std::max(reducedEnergy,energyLowLimit);
118 
119  // Helium ion case
120  if( Zi < 2.5 ) {
121 
122  static const G4double c[6] = {0.2865,0.1266,-0.001429,0.02402,-0.01135,0.001475};
123 
124  G4double Q = std::max(0.0,G4Log(reducedEnergy*massFactor));
125  G4double x = c[0];
126  G4double y = 1.0;
127  for (G4int i=1; i<6; i++) {
128  y *= Q;
129  x += y * c[i] ;
130  }
131  G4double ex;
132  if(x < 0.2) ex = x * (1 - 0.5*x);
133  else ex = 1. - G4Exp(-x);
134 
135  G4double tq = 7.6 - Q;
136  G4double tq2= tq*tq;
137  G4double tt = ( 0.007 + 0.00005 * z );
138  if(tq2 < 0.2) tt *= (1.0 - tq2 + 0.5*tq2*tq2);
139  else tt *= G4Exp(-tq2);
140 
141  effCharge = charge*(1.0 + tt) * std::sqrt(ex);
142 
143  // Heavy ion case
144  } else {
145 
146  G4double y;
147  // = g4pow->A13(Zi);
148  //G4double z23 = y*y;
149  G4double zi13 = g4pow->A13(Zi);
150  G4double zi23 = zi13*zi13;
151  // G4double e = std::max(reducedEnergy,energyBohr/z23);
152  //G4double e = reducedEnergy;
153 
154  // v1 is ion velocity in vF unit
155  G4double eF = material->GetIonisation()->GetFermiEnergy();
156  G4double v1sq = reducedEnergy/eF;
157  G4double vFsq = eF/energyBohr;
158  G4double vF = std::sqrt(eF/energyBohr);
159 
160  // Faster than Fermi velocity
161  if ( v1sq > 1.0 ) {
162  y = vF * std::sqrt(v1sq) * ( 1.0 + 0.2/v1sq ) / zi23 ;
163 
164  // Slower than Fermi velocity
165  } else {
166  y = 0.692308 * vF * (1.0 + 0.666666*v1sq + v1sq*v1sq/15.0) / zi23 ;
167  }
168 
169  G4double q;
170  G4double y3 = std::pow(y, 0.3) ;
171  // G4cout<<"y= "<<y<<" y3= "<<y3<<" v1= "<<v1<<" vF= "<<vF<<G4endl;
172  q = 1.0 - G4Exp( 0.803*y3 - 1.3167*y3*y3 - 0.38157*y - 0.008983*y*y ) ;
173 
174  //y *= 0.77;
175  //y *= (0.75 + 0.52/Zi);
176 
177  //if( y < 0.2 ) q = y*(1.0 - 0.5*y);
178  //else q = 1.0 - G4Exp(-y);
179 
180  G4double qmin = minCharge/Zi;
181  if(q < qmin) q = qmin;
182 
183  effCharge = q*charge;
184 
185  /*
186  G4double x1 = 1.0*effCharge*(1.0 - 0.132*G4Log(y))/(y*std::sqrt(z));
187  G4double x2 = 0.1*effCharge*effCharge*energyBohr/reducedEnergy;
188 
189  chargeCorrection = 1.0 + x1 - x2;
190 
191  G4cout << "x1= "<<x1<<" x2= "<< x2<<" corr= "<<chargeCorrection<<G4endl;
192  */
193 
194  G4double tq = 7.6 - G4Log(reducedEnergy/keV);
195  G4double tq2= tq*tq;
196  G4double sq = 1.0 + ( 0.18 + 0.0015 * z )*G4Exp(-tq2)/ (Zi*Zi);
197 
198  // G4cout << "sq= " << sq << G4endl;
199 
200  // Screen length according to
201  // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
202  // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
203 
204  G4double lambda = 10.0 * vF *g4pow->A23(1.0 - q)/ (zi13 * (6.0 + q));
205 
206  G4double lambda2 = lambda*lambda;
207 
208  G4double xx = (0.5/q - 0.5)*G4Log(1.0 + lambda2)/vFsq;
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 
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double z
Definition: TRTMaterials.hh:39
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
G4double GetZeffective() const
G4double GetFermiEnergy() const
string material
Definition: eplot.py:19
G4double A23(G4double A) const
Definition: G4Pow.hh:159
float proton_mass_c2
Definition: hepunit.py:275
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4double GetPDGMass() const
G4double A13(G4double A) const
Definition: G4Pow.hh:134
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76
float amu_c2
Definition: hepunit.py:277
G4double GetPDGCharge() const