Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4ionEffectiveCharge
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 07.05.2002
37//
38// Modifications:
39// 12.09.2004 Set low energy limit to 1 keV (V.Ivanchenko)
40// 25.01.2005 Add protection - min Charge 0.1 eplus (V.Ivanchenko)
41// 28.04.2006 Set upper energy limit to 50 MeV (V.Ivanchenko)
42// 23.05.2006 Set upper energy limit to Z*10 MeV (V.Ivanchenko)
43// 15.08.2006 Add protection for not defined material (V.Ivanchenko)
44// 27-09-2007 Use Fermi energy from material, optimazed formulas (V.Ivanchenko)
45//
46
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
54#include "G4SystemOfUnits.hh"
55#include "G4UnitsTable.hh"
56#include "G4Material.hh"
57#include "G4NistManager.hh"
58#include "G4Log.hh"
59#include "G4Exp.hh"
60#include "G4Pow.hh"
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65{
66 chargeCorrection = 1.0;
67 energyHighLimit = 20.0*CLHEP::MeV;
68 energyLowLimit = 1.0*CLHEP::keV;
69 energyBohr = 25.*CLHEP::keV;
70 massFactor = CLHEP::amu_c2/(CLHEP::proton_mass_c2*CLHEP::keV);
71 minCharge = 1.0;
72 lastKinEnergy = 0.0;
73 effCharge = CLHEP::eplus;
74 inveplus = 1.0/CLHEP::eplus;
75 g4calc = G4Pow::GetInstance();
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81 const G4Material* material,
82 G4double kineticEnergy)
83{
84 if(p == lastPart && material == lastMat && kineticEnergy == lastKinEnergy)
85 return effCharge;
86
87 lastPart = p;
88 lastMat = material;
89 lastKinEnergy = kineticEnergy;
90
91 G4double mass = p->GetPDGMass();
92 effCharge = p->GetPDGCharge();
93 G4int Zi = G4lrint(effCharge*inveplus);
94 chargeCorrection = 1.0;
95 if(Zi <= 1) { return effCharge; }
96
97 // The aproximation of ion effective charge from:
98 // J.F.Ziegler, J.P. Biersack, U. Littmark
99 // The Stopping and Range of Ions in Matter,
100 // Vol.1, Pergamon Press, 1985
101 // Fast ions or hadrons
102 G4double reducedEnergy = kineticEnergy * CLHEP::proton_mass_c2/mass;
103
104 //G4cout << "e= " << reducedEnergy << " Zi= " << Zi << " "
105 //<< material->GetName() << G4endl;
106
107 if(reducedEnergy > effCharge*energyHighLimit ) {
108 return effCharge;
109 }
110 G4double z = material->GetIonisation()->GetZeffective();
111 reducedEnergy = std::max(reducedEnergy,energyLowLimit);
112
113 // Helium ion case
114 if( Zi <= 2 ) {
115
116 static const G4double c[6] =
117 {0.2865,0.1266,-0.001429,0.02402,-0.01135,0.001475};
118
119 G4double Q = std::max(0.0,G4Log(reducedEnergy*massFactor));
120 G4double x = c[0];
121 G4double y = 1.0;
122 for (G4int i=1; i<6; ++i) {
123 y *= Q;
124 x += y * c[i] ;
125 }
126 G4double ex = (x < 0.2) ? x * (1 - 0.5*x) : 1. - G4Exp(-x);
127
128 G4double tq = 7.6 - Q;
129 G4double tq2= tq*tq;
130 G4double tt = ( 0.007 + 0.00005 * z );
131 if(tq2 < 0.2) { tt *= (1.0 - tq2 + 0.5*tq2*tq2); }
132 else { tt *= G4Exp(-tq2); }
133
134 effCharge *= (1.0 + tt) * std::sqrt(ex);
135
136 // Heavy ion case
137 } else {
138
139 G4double zi13 = g4calc->Z13(Zi);
140 G4double zi23 = zi13*zi13;
141
142 // v1 is ion velocity in vF unit
143 G4double eF = material->GetIonisation()->GetFermiEnergy();
144 G4double v1sq = reducedEnergy/eF;
145 G4double vFsq = eF/energyBohr;
146 G4double vF = std::sqrt(eF/energyBohr);
147
148 G4double y = ( v1sq > 1.0 )
149 // Faster than Fermi velocity
150 ? vF * std::sqrt(v1sq) * ( 1.0 + 0.2/v1sq ) / zi23
151 // Slower than Fermi velocity
152 : 0.692308 * vF * (1.0 + 0.666666*v1sq + v1sq*v1sq/15.0) / zi23;
153
154 G4double y3 = G4Exp(0.3*G4Log(y));
155 // G4cout<<"y= "<<y<<" y3= "<<y3<<" v1= "<<v1<<" vF= "<<vF<<G4endl;
156 G4double q = std::max(1.0 - G4Exp( 0.803*y3 - 1.3167*y3*y3 - 0.38157*y
157 - 0.008983*y*y), minCharge/effCharge);
158
159 // compute charge correction
160 G4double tq = 7.6 - G4Log(reducedEnergy/CLHEP::keV);
161 G4double tq2= tq*tq;
162 G4double sq = 1.0 + ( 0.18 + 0.0015 * z )*G4Exp(-tq2)/ (Zi*Zi);
163 // G4cout << "sq= " << sq << G4endl;
164
165 // Screen length according to
166 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
167 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
168
169 G4double lambda = 10.0 * vF *g4calc->A23(1.0 - q)/ (zi13 * (6.0 + q));
170 G4double lambda2 = lambda*lambda;
171 G4double xx = (0.5/q - 0.5)*G4Log(1.0 + lambda2)/vFsq;
172
173 effCharge *= q;
174 chargeCorrection = sq * (1.0 + xx);
175 }
176 // G4cout << "G4ionEffectiveCharge: charge= " << charge << " q= " << q
177 // << " chargeCor= " << chargeCorrection
178 // << " e(MeV)= " << kineticEnergy/MeV << G4endl;
179 return effCharge;
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double GetZeffective() const
G4double GetFermiEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4double A23(G4double A) const
Definition: G4Pow.hh:131
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
int G4lrint(double ad)
Definition: templates.hh:134