Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ionEffectiveCharge.cc
이 파일의 문서화 페이지로 가기
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 108386 2018-02-09 15:38:32Z 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;
72  minCharge = 1.0;
73  lastPart = 0;
74  lastMat = 0;
75  lastKinEnergy = 0.0;
76  effCharge = eplus;
77  inveplus = 1.0/CLHEP::eplus;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 
84 {}
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87 
89  const G4Material* material,
91 {
92  if(p == lastPart && material == lastMat && kineticEnergy == lastKinEnergy)
93  return effCharge;
94 
95  lastPart = p;
96  lastMat = material;
98 
99  G4double mass = p->GetPDGMass();
100  G4double charge = p->GetPDGCharge();
101  G4double Zi = charge*inveplus;
102 
103  chargeCorrection = 1.0;
104  effCharge = charge;
105 
106  // The aproximation of ion effective charge from:
107  // J.F.Ziegler, J.P. Biersack, U. Littmark
108  // The Stopping and Range of Ions in Matter,
109  // Vol.1, Pergamon Press, 1985
110  // Fast ions or hadrons
111  G4double reducedEnergy = kineticEnergy * proton_mass_c2/mass ;
112 
113  //G4cout << "e= " << reducedEnergy << " Zi= " << Zi << " "
114  //<< material->GetName() << G4endl;
115 
116  if(Zi < 1.5 || !material || reducedEnergy > Zi*energyHighLimit ) {
117  return charge;
118  }
119  G4double z = material->GetIonisation()->GetZeffective();
120  reducedEnergy = std::max(reducedEnergy,energyLowLimit);
121 
122  // Helium ion case
123  if( Zi < 2.5 ) {
124 
125  static const G4double c[6] =
126  {0.2865,0.1266,-0.001429,0.02402,-0.01135,0.001475};
127 
128  G4double Q = std::max(0.0,G4Log(reducedEnergy*massFactor));
129  G4double x = c[0];
130  G4double y = 1.0;
131  for (G4int i=1; i<6; ++i) {
132  y *= Q;
133  x += y * c[i] ;
134  }
135  G4double ex;
136  if(x < 0.2) { ex = x * (1 - 0.5*x); }
137  else { ex = 1. - G4Exp(-x); }
138 
139  G4double tq = 7.6 - Q;
140  G4double tq2= tq*tq;
141  G4double tt = ( 0.007 + 0.00005 * z );
142  if(tq2 < 0.2) { tt *= (1.0 - tq2 + 0.5*tq2*tq2); }
143  else { tt *= G4Exp(-tq2); }
144 
145  effCharge = charge*(1.0 + tt) * std::sqrt(ex);
146 
147  // Heavy ion case
148  } else {
149 
150  G4double y;
151  G4double zi13 = g4calc->A13(Zi);
152  G4double zi23 = zi13*zi13;
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  q = std::max(q, minCharge/Zi);
174 
175  effCharge = q*charge;
176 
177  G4double tq = 7.6 - G4Log(reducedEnergy/keV);
178  G4double tq2= tq*tq;
179  G4double sq = 1.0 + ( 0.18 + 0.0015 * z )*G4Exp(-tq2)/ (Zi*Zi);
180 
181  // G4cout << "sq= " << sq << G4endl;
182 
183  // Screen length according to
184  // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
185  // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
186 
187  G4double lambda = 10.0 * vF *g4calc->A23(1.0 - q)/ (zi13 * (6.0 + q));
188 
189  G4double lambda2 = lambda*lambda;
190 
191  G4double xx = (0.5/q - 0.5)*G4Log(1.0 + lambda2)/vFsq;
192 
193  chargeCorrection = sq * (1.0 + xx);
194 
195  }
196  // G4cout << "G4ionEffectiveCharge: charge= " << charge << " q= " << q
197  // << " chargeCor= " << chargeCorrection
198  // << " e(MeV)= " << kineticEnergy/MeV << G4endl;
199  return effCharge;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203 
204 
Float_t x
Definition: compare.C:6
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Double_t xx
G4double A13(G4double A) const
Definition: G4Pow.cc:138
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double keV
Definition: G4SIunits.hh:216
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
Double_t z
G4double A23(G4double A) const
Definition: G4Pow.hh:143
const G4Material * lastMat
G4double GetPDGCharge() const
G4double GetPDGMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double proton_mass_c2
double G4double
Definition: G4Types.hh:76
static constexpr double amu_c2
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
G4double GetFermiEnergy() const
static constexpr double eplus
Definition: G4SIunits.hh:199
int G4int
Definition: G4Types.hh:78
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:227
G4double GetZeffective() const
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
static double Q[]
const G4ParticleDefinition * lastPart
static constexpr double eplus