Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4hIonEffChargeSquare.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 //
27 // -------------------------------------------------------------------
28 //
29 // GEANT4 Class file
30 //
31 //
32 // File name: G4hIonEffChargeSquare
33 //
34 // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
35 //
36 // Creation date: 20 July 2000
37 //
38 // Modifications:
39 // 20/07/2000 V.Ivanchenko First implementation
40 // 18/06/2001 V.Ivanchenko Continuation for eff.charge (small change of y)
41 // 08/10/2002 V.Ivanchenko The charge of the nucleus is used not charge of
42 // DynamicParticle
43 //
44 // Class Description:
45 //
46 // Ion effective charge model
47 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
48 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
49 //
50 // Class Description: End
51 //
52 // -------------------------------------------------------------------
53 //
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 
56 #include "G4hIonEffChargeSquare.hh"
57 #include "G4PhysicalConstants.hh"
58 #include "G4SystemOfUnits.hh"
59 #include "G4DynamicParticle.hh"
60 #include "G4ParticleDefinition.hh"
61 #include "G4Material.hh"
62 #include "G4Element.hh"
63 #include "G4Exp.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
68  : G4VLowEnergyModel(name),
69  theHeMassAMU(4.0026)
70 {;}
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73 
75 {;}
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80  const G4Material* material)
81 {
82  G4double energy = particle->GetKineticEnergy() ;
83  G4double particleMass = particle->GetMass() ;
84  G4double charge = (particle->GetDefinition()->GetPDGCharge())/eplus ;
85 
86  G4double q = IonEffChargeSquare(material,energy,particleMass,charge) ;
87 
88  return q ;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94  const G4Material* material,
96 {
97  // SetRateMass(aParticle) ;
98  G4double particleMass = aParticle->GetPDGMass() ;
99  G4double charge = (aParticle->GetPDGCharge())/eplus ;
100 
101  G4double q = IonEffChargeSquare(material,kineticEnergy,particleMass,charge) ;
102 
103  return q ;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109  const G4ParticleDefinition* ,
110  const G4Material* ) const
111 {
112  return 1.0*TeV ;
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
118  const G4ParticleDefinition* ,
119  const G4Material* ) const
120 {
121  return 0.0 ;
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125 
127  const G4ParticleDefinition* ) const
128 {
129  return 1.0*TeV ;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
135  const G4ParticleDefinition* ) const
136 {
137  return 0.0 ;
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
143  const G4Material* ) const
144 {
145  return true ;
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149 
151  const G4Material* ) const
152 {
153  return true ;
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157 
159  const G4Material* material,
161  G4double particleMass,
162  G4double ionCharge) const
163 {
164  // The aproximation of ion effective charge from:
165  // J.F.Ziegler, J.P. Biersack, U. Littmark
166  // The Stopping and Range of Ions in Matter,
167  // Vol.1, Pergamon Press, 1985
168 
169  // Fast ions or hadrons
170  G4double reducedEnergy = kineticEnergy * proton_mass_c2/particleMass ;
171  if(reducedEnergy < 1.0*keV) reducedEnergy = 1.0*keV;
172  if( (reducedEnergy > ionCharge * 10.0 * MeV) ||
173  (ionCharge < 1.5) ) return ionCharge*ionCharge ;
174 
175  static const G4double vFermi[92] = {
176  1.0309, 0.15976, 0.59782, 1.0781, 1.0486, 1.0, 1.058, 0.93942, 0.74562, 0.3424,
177  0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712,
178  0.81707, 0.9943, 1.1423, 1.2381, 1.1222, 0.92705, 1.0047, 1.2, 1.0661, 0.97411,
179  0.84912, 0.95, 1.0903, 1.0429, 0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207,
180  1.029, 1.2542, 1.122, 1.1241, 1.0882, 1.2709, 1.2542, 0.90094, 0.74093, 0.86054,
181  0.93155, 1.0047, 0.55379, 0.43289, 0.32636, 0.5131, 0.695, 0.72591, 0.71202, 0.67413,
182  0.71418, 0.71453, 0.5911, 0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884,
183  0.71801, 0.83048, 1.1222, 1.2381, 1.045, 1.0733, 1.0953, 1.2381, 1.2879, 0.78654,
184  0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065,
185  1.9578, 1.0257} ;
186 
187  static const G4double lFactor[92] = {
188  1.0, 1.0, 1.1, 1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9,
189  0.82, 0.81, 0.83, 0.88, 1.0, 0.95, 0.97, 0.99, 0.98, 0.97,
190  0.98, 0.97, 0.96, 0.93, 0.91, 0.9, 0.88, 0.9, 0.9, 0.9,
191  0.9, 0.85, 0.9, 0.9, 0.91, 0.92, 0.9, 0.9, 0.9, 0.9,
192  0.9, 0.88, 0.9, 0.88, 0.88, 0.9, 0.9, 0.88, 0.9, 0.9,
193  0.9, 0.9, 0.96, 1.2, 0.9, 0.88, 0.88, 0.85, 0.9, 0.9,
194  0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1, 1.08, 1.08,
195  1.08, 1.08, 1.09, 1.09, 1.1, 1.11, 1.12, 1.13, 1.14, 1.15,
196  1.17, 1.2, 1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16,
197  1.16, 1.16} ;
198 
199  static const G4double c[6] = {0.2865, 0.1266, -0.001429,
200  0.02402,-0.01135, 0.001475} ;
201 
202  // get elements in the actual material,
203  const G4ElementVector* theElementVector = material->GetElementVector() ;
204  const G4double* theAtomicNumDensityVector =
205  material->GetAtomicNumDensityVector() ;
206  const G4int NumberOfElements = material->GetNumberOfElements() ;
207 
208  // loop for the elements in the material
209  // to find out average values Z, vF, lF
210  G4double z = 0.0, vF = 0.0, lF = 0.0, norm = 0.0 ;
211 
212  if( 1 == NumberOfElements ) {
213  z = material->GetZ() ;
214  G4int iz = G4int(z) - 1 ;
215  if(iz < 0) iz = 0 ;
216  else if(iz > 91) iz = 91 ;
217  vF = vFermi[iz] ;
218  lF = lFactor[iz] ;
219 
220  } else {
221  for (G4int iel=0; iel<NumberOfElements; iel++)
222  {
223  const G4Element* element = (*theElementVector)[iel] ;
224  G4double z2 = element->GetZ() ;
225  const G4double weight = theAtomicNumDensityVector[iel] ;
226  norm += weight ;
227  z += z2 * weight ;
228  G4int iz = G4int(z2) - 1 ;
229  if(iz < 0) iz = 0 ;
230  else if(iz > 91) iz =91 ;
231  vF += vFermi[iz] * weight ;
232  lF += lFactor[iz] * weight ;
233  }
234  z /= norm ;
235  vF /= norm ;
236  lF /= norm ;
237  }
238 
239  // Helium ion case
240  if( ionCharge < 2.5 ) {
241 
242  G4double e = std::log(std::max(1.0, kineticEnergy / (keV*theHeMassAMU) )) ;
243  G4double x = c[0] ;
244  G4double y = 1.0 ;
245  for (G4int i=1; i<6; i++) {
246  y *= e ;
247  x += y * c[i] ;
248  }
249  G4double q = 7.6 - e ;
250  q = 1.0 + ( 0.007 + 0.00005 * z ) * G4Exp( -q*q ) ;
251  return 4.0 * q * q * (1.0 - G4Exp(-x)) ;
252 
253  // Heavy ion case
254  } else {
255 
256  // v1 is ion velocity in vF unit
257  G4double v1 = std::sqrt( reducedEnergy / (25.0 * keV) )/ vF ;
258  G4double y ;
259  G4double z13 = std::pow(ionCharge, 0.3333) ;
260 
261  // Faster than Fermi velocity
262  if ( v1 > 1.0 ) {
263  y = vF * v1 * ( 1.0 + 0.2 / (v1*v1) ) / (z13*z13) ;
264 
265  // Slower than Fermi velocity
266  } else {
267  y = 0.6923 * vF * (1.0 + 2.0*v1*v1/3.0 + v1*v1*v1*v1/15.0) / (z13*z13) ;
268  }
269 
270  G4double y3 = std::pow(y, 0.3) ;
271  G4double q = 1.0 - G4Exp( 0.803*y3 - 1.3167*y3*y3 -
272  0.38157*y - 0.008983*y*y ) ;
273  if( q < 0.0 ) q = 0.0 ;
274 
275  G4double sLocal = 7.6 - std::log(std::max(1.0, reducedEnergy/keV)) ;
276  sLocal = 1.0 + ( 0.18 + 0.0015 * z ) * G4Exp( -sLocal*sLocal )/ (ionCharge*ionCharge) ;
277 
278  // Screen length according to
279  // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
280  // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
281 
282  G4double lambda = 10.0 * vF * std::pow(1.0-q, 0.6667) / (z13 * (6.0 + q)) ;
283  G4double qeff = ionCharge * sLocal *
284  ( q + 0.5*(1.0-q) * std::log(1.0 + lambda*lambda) / (vF*vF) ) ;
285  if( 0.1 > qeff ) qeff = 0.1 ;
286  return qeff*qeff ;
287  }
288 }
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
const XML_Char * name
Definition: expat.h:151
static constexpr double MeV
Definition: G4SIunits.hh:214
G4hIonEffChargeSquare(const G4String &name)
static constexpr double keV
Definition: G4SIunits.hh:216
Float_t y
Definition: compare.C:6
G4double HighEnergyLimit(const G4ParticleDefinition *aParticle, const G4Material *material) const
Double_t z
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
G4double GetPDGCharge() const
G4double GetPDGMass() const
static constexpr double proton_mass_c2
static constexpr double TeV
Definition: G4SIunits.hh:218
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
double energy
Definition: plottest35.C:25
G4double GetZ() const
Definition: G4Material.cc:629
G4double LowEnergyLimit(const G4ParticleDefinition *aParticle, const G4Material *material) const
static constexpr double eplus
Definition: G4SIunits.hh:199
double weight
Definition: plottest35.C:25
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
Float_t norm
G4double IonEffChargeSquare(const G4Material *material, G4double kineticEnergy, G4double particleMass, G4double ionCharge) const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetMass() const
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:217
G4double TheValue(const G4DynamicParticle *particle, const G4Material *material)
G4bool IsInCharge(const G4DynamicParticle *particle, const G4Material *material) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:187