Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ICRU49NuclearStoppingModel.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: G4ICRU49NuclearStoppingModel.cc 100399 2016-10-20 07:38:12Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4ICRU49NuclearStoppingModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 09.04.2008 from G4MuMscModel
38 //
39 // Modifications:
40 //
41 //
42 // Class Description:
43 //
44 // Implementation of the model of ICRU'49 nuclear stopping
45 
46 // -------------------------------------------------------------------
47 //
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
53 #include "G4PhysicalConstants.hh"
54 #include "G4SystemOfUnits.hh"
55 #include "Randomize.hh"
56 #include "G4LossTableManager.hh"
58 #include "G4ElementVector.hh"
59 #include "G4ProductionCutsTable.hh"
60 #include "G4Step.hh"
61 #include "G4Pow.hh"
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
65 using namespace std;
66 
68  : G4VEmModel(nam),lossFlucFlag(false)
69 {
70  theZieglerFactor = eV*cm2*1.0e-15;
71  g4calc = G4Pow::GetInstance();
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75 
77 {}
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82  const G4DataVector&)
83 {}
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
87 void
89  std::vector<G4DynamicParticle*>*,
90  const G4MaterialCutsCouple*,
91  const G4DynamicParticle*,
93 {}
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
97 G4double
99  const G4Material* mat,
100  const G4ParticleDefinition* p,
101  G4double kinEnergy,
102  G4double)
103 {
104  G4double nloss = 0.0;
105  if(kinEnergy <= 0.0) { return nloss; }
106 
107  // projectile
108  G4double mass1 = p->GetPDGMass();
109  G4double z1 = std::fabs(p->GetPDGCharge()/eplus);
110 
111  if(kinEnergy*proton_mass_c2/mass1 > z1*z1*MeV) { return nloss; }
112 
113  // Projectile nucleus
114  mass1 /= amu_c2;
115 
116  // loop for the elements in the material
117  G4int numberOfElements = mat->GetNumberOfElements();
118  const G4ElementVector* theElementVector = mat->GetElementVector();
119  const G4double* atomDensity = mat->GetAtomicNumDensityVector();
120 
121  for (G4int iel=0; iel<numberOfElements; iel++) {
122  const G4Element* element = (*theElementVector)[iel] ;
123  G4double z2 = element->GetZ();
124  G4double mass2 = element->GetN();
125  nloss += (NuclearStoppingPower(kinEnergy, z1, z2, mass1, mass2))
126  * atomDensity[iel] ;
127  }
128  nloss *= theZieglerFactor;
129  return nloss;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
134 G4double
135 G4ICRU49NuclearStoppingModel::NuclearStoppingPower(G4double kineticEnergy,
136  G4double z1, G4double z2,
137  G4double mass1, G4double mass2)
138 {
139  G4double energy = kineticEnergy/keV ; // energy in keV
140  G4double nloss = 0.0;
141  G4double z12 = z1*z2;
142  G4int iz1 = G4lrint(z1);
143  G4int iz2 = G4lrint(z2);
144 
145  G4double rm;
146  if(iz1 > 1) { rm = (mass1 + mass2)*(g4calc->Z23(iz1) + g4calc->Z23(iz2)); }
147  else { rm = (mass1 + mass2)*g4calc->Z13(iz2); }
148 
149  G4double er = 32.536 * mass2 * energy / ( z12 * rm ) ; // reduced energy
150 
151  static const G4double nuca[104][2] = {
152  { 1.0E+8, 5.831E-8},
153  { 8.0E+7, 7.288E-8},
154  { 6.0E+7, 9.719E-8},
155  { 5.0E+7, 1.166E-7},
156  { 4.0E+7, 1.457E-7},
157  { 3.0E+7, 1.942E-7},
158  { 2.0E+7, 2.916E-7},
159  { 1.5E+7, 3.887E-7},
160 
161  { 1.0E+7, 5.833E-7},
162  { 8.0E+6, 7.287E-7},
163  { 6.0E+6, 9.712E-7},
164  { 5.0E+6, 1.166E-6},
165  { 4.0E+6, 1.457E-6},
166  { 3.0E+6, 1.941E-6},
167  { 2.0E+6, 2.911E-6},
168  { 1.5E+6, 3.878E-6},
169 
170  { 1.0E+6, 5.810E-6},
171  { 8.0E+5, 7.262E-6},
172  { 6.0E+5, 9.663E-6},
173  { 5.0E+5, 1.157E-5},
174  { 4.0E+5, 1.442E-5},
175  { 3.0E+5, 1.913E-5},
176  { 2.0E+5, 2.845E-5},
177  { 1.5E+5, 3.762E-5},
178 
179  { 1.0E+5, 5.554E-5},
180  { 8.0E+4, 6.866E-5},
181  { 6.0E+4, 9.020E-5},
182  { 5.0E+4, 1.070E-4},
183  { 4.0E+4, 1.319E-4},
184  { 3.0E+4, 1.722E-4},
185  { 2.0E+4, 2.499E-4},
186  { 1.5E+4, 3.248E-4},
187 
188  { 1.0E+4, 4.688E-4},
189  { 8.0E+3, 5.729E-4},
190  { 6.0E+3, 7.411E-4},
191  { 5.0E+3, 8.718E-4},
192  { 4.0E+3, 1.063E-3},
193  { 3.0E+3, 1.370E-3},
194  { 2.0E+3, 1.955E-3},
195  { 1.5E+3, 2.511E-3},
196 
197  { 1.0E+3, 3.563E-3},
198  { 8.0E+2, 4.314E-3},
199  { 6.0E+2, 5.511E-3},
200  { 5.0E+2, 6.430E-3},
201  { 4.0E+2, 7.756E-3},
202  { 3.0E+2, 9.855E-3},
203  { 2.0E+2, 1.375E-2},
204  { 1.5E+2, 1.736E-2},
205 
206  { 1.0E+2, 2.395E-2},
207  { 8.0E+1, 2.850E-2},
208  { 6.0E+1, 3.552E-2},
209  { 5.0E+1, 4.073E-2},
210  { 4.0E+1, 4.802E-2},
211  { 3.0E+1, 5.904E-2},
212  { 1.5E+1, 9.426E-2},
213 
214  { 1.0E+1, 1.210E-1},
215  { 8.0E+0, 1.377E-1},
216  { 6.0E+0, 1.611E-1},
217  { 5.0E+0, 1.768E-1},
218  { 4.0E+0, 1.968E-1},
219  { 3.0E+0, 2.235E-1},
220  { 2.0E+0, 2.613E-1},
221  { 1.5E+0, 2.871E-1},
222 
223  { 1.0E+0, 3.199E-1},
224  { 8.0E-1, 3.354E-1},
225  { 6.0E-1, 3.523E-1},
226  { 5.0E-1, 3.609E-1},
227  { 4.0E-1, 3.693E-1},
228  { 3.0E-1, 3.766E-1},
229  { 2.0E-1, 3.803E-1},
230  { 1.5E-1, 3.788E-1},
231 
232  { 1.0E-1, 3.711E-1},
233  { 8.0E-2, 3.644E-1},
234  { 6.0E-2, 3.530E-1},
235  { 5.0E-2, 3.444E-1},
236  { 4.0E-2, 3.323E-1},
237  { 3.0E-2, 3.144E-1},
238  { 2.0E-2, 2.854E-1},
239  { 1.5E-2, 2.629E-1},
240 
241  { 1.0E-2, 2.298E-1},
242  { 8.0E-3, 2.115E-1},
243  { 6.0E-3, 1.883E-1},
244  { 5.0E-3, 1.741E-1},
245  { 4.0E-3, 1.574E-1},
246  { 3.0E-3, 1.372E-1},
247  { 2.0E-3, 1.116E-1},
248  { 1.5E-3, 9.559E-2},
249 
250  { 1.0E-3, 7.601E-2},
251  { 8.0E-4, 6.668E-2},
252  { 6.0E-4, 5.605E-2},
253  { 5.0E-4, 5.008E-2},
254  { 4.0E-4, 4.352E-2},
255  { 3.0E-4, 3.617E-2},
256  { 2.0E-4, 2.768E-2},
257  { 1.5E-4, 2.279E-2},
258 
259  { 1.0E-4, 1.723E-2},
260  { 8.0E-5, 1.473E-2},
261  { 6.0E-5, 1.200E-2},
262  { 5.0E-5, 1.052E-2},
263  { 4.0E-5, 8.950E-3},
264  { 3.0E-5, 7.246E-3},
265  { 2.0E-5, 5.358E-3},
266  { 1.5E-5, 4.313E-3},
267  { 0.0, 3.166E-3}
268  };
269 
270  if (er >= nuca[0][0]) { nloss = nuca[0][1]; }
271  else {
272  // the table is inverse in energy
273  for (G4int i=102; i>=0; --i) {
274  G4double edi = nuca[i][0];
275  if (er <= edi) {
276  G4double edi1 = nuca[i+1][0];
277  G4double ai = nuca[i][1];
278  G4double ai1 = nuca[i+1][1];
279  nloss = (ai - ai1)*(er - edi1)/(edi - edi1) + ai1;
280  break;
281  }
282  }
283  }
284 
285  // Stragling
286  if(lossFlucFlag) {
287  G4double sig = 4.0 * mass1 * mass2 / ((mass1 + mass2)*(mass1 + mass2)*
288  (4.0 + 0.197/(er*er) + 6.584/er));
289 
290  nloss *= G4RandGauss::shoot(1.0,sig);
291  lossFlucFlag = false;
292  }
293 
294  nloss *= 8.462 * z12 * mass1 / rm; // Return to [ev/(10^15 atoms/cm^2]
295 
296  nloss = std::max(nloss, 0.0);
297  return nloss;
298 }
299 
300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
ThreeVector shoot(const G4int Ap, const G4int Af)
std::vector< G4Element * > G4ElementVector
static constexpr double proton_mass_c2
static constexpr double cm2
Definition: G4SIunits.hh:120
G4double GetN() const
Definition: G4Element.hh:135
const char * p
Definition: xmltok.h:285
G4double GetZ() const
Definition: G4Element.hh:131
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
static constexpr double eplus
Definition: G4SIunits.hh:199
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double amu_c2
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double) final
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX) final
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
G4double Z23(G4int Z) const
Definition: G4Pow.hh:154
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4ICRU49NuclearStoppingModel(const G4String &nam="ICRU49NucStopping")
static constexpr double keV
Definition: G4SIunits.hh:216