Geant4-11
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;
71 minCharge = 1.0;
72 lastKinEnergy = 0.0;
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;
89 lastKinEnergy = kineticEnergy;
90
91 G4double mass = p->GetPDGMass();
94 chargeCorrection = 1.0;
95
96 // The aproximation of ion effective charge from:
97 // J.F.Ziegler, J.P. Biersack, U. Littmark
98 // The Stopping and Range of Ions in Matter,
99 // Vol.1, Pergamon Press, 1985
100 // Fast ions or hadrons
101 G4double reducedEnergy = kineticEnergy * CLHEP::proton_mass_c2/mass;
102
103 //G4cout << "e= " << reducedEnergy << " Zi= " << Zi << " "
104 //<< material->GetName() << G4endl;
105
106 if(Zi <= 1 || reducedEnergy > effCharge*energyHighLimit ) {
107 return effCharge;
108 }
109 G4double z = material->GetIonisation()->GetZeffective();
110 reducedEnergy = std::max(reducedEnergy,energyLowLimit);
111
112 // Helium ion case
113 if( Zi <= 2 ) {
114
115 static const G4double c[6] =
116 {0.2865,0.1266,-0.001429,0.02402,-0.01135,0.001475};
117
118 G4double Q = std::max(0.0,G4Log(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 = (x < 0.2) ? x * (1 - 0.5*x) : 1. - G4Exp(-x);
126
127 G4double tq = 7.6 - Q;
128 G4double tq2= tq*tq;
129 G4double tt = ( 0.007 + 0.00005 * z );
130 if(tq2 < 0.2) { tt *= (1.0 - tq2 + 0.5*tq2*tq2); }
131 else { tt *= G4Exp(-tq2); }
132
133 effCharge *= (1.0 + tt) * std::sqrt(ex);
134
135 // Heavy ion case
136 } else {
137
138 G4double zi13 = g4calc->Z13(Zi);
139 G4double zi23 = zi13*zi13;
140
141 // v1 is ion velocity in vF unit
142 G4double eF = material->GetIonisation()->GetFermiEnergy();
143 G4double v1sq = reducedEnergy/eF;
144 G4double vFsq = eF/energyBohr;
145 G4double vF = std::sqrt(eF/energyBohr);
146
147 G4double y = ( v1sq > 1.0 )
148 // Faster than Fermi velocity
149 ? vF * std::sqrt(v1sq) * ( 1.0 + 0.2/v1sq ) / zi23
150 // Slower than Fermi velocity
151 : 0.692308 * vF * (1.0 + 0.666666*v1sq + v1sq*v1sq/15.0) / zi23;
152
153 G4double y3 = G4Exp(0.3*G4Log(y));
154 // G4cout<<"y= "<<y<<" y3= "<<y3<<" v1= "<<v1<<" vF= "<<vF<<G4endl;
155 G4double q = std::max(1.0 - G4Exp( 0.803*y3 - 1.3167*y3*y3 - 0.38157*y
156 - 0.008983*y*y), minCharge/effCharge);
157
158 // compute charge correction
159 G4double tq = 7.6 - G4Log(reducedEnergy/CLHEP::keV);
160 G4double tq2= tq*tq;
161 G4double sq = 1.0 + ( 0.18 + 0.0015 * z )*G4Exp(-tq2)/ (Zi*Zi);
162 // G4cout << "sq= " << sq << G4endl;
163
164 // Screen length according to
165 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
166 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
167
168 G4double lambda = 10.0 * vF *g4calc->A23(1.0 - q)/ (zi13 * (6.0 + q));
169 G4double lambda2 = lambda*lambda;
170 G4double xx = (0.5/q - 0.5)*G4Log(1.0 + lambda2)/vFsq;
171
172 effCharge *= q;
173 chargeCorrection = sq * (1.0 + xx);
174 }
175 // G4cout << "G4ionEffectiveCharge: charge= " << charge << " q= " << q
176 // << " chargeCor= " << chargeCorrection
177 // << " e(MeV)= " << kineticEnergy/MeV << G4endl;
178 return effCharge;
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
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)
const G4Material * lastMat
const G4ParticleDefinition * lastPart
static constexpr double eplus
static constexpr double amu_c2
static constexpr double proton_mass_c2
static constexpr double keV
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
string material
Definition: eplot.py:19
static double Q[]
int G4lrint(double ad)
Definition: templates.hh:134