Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
54#include "G4SystemOfUnits.hh"
55#include "Randomize.hh"
56#include "G4LossTableManager.hh"
58#include "G4ElementVector.hh"
60#include "G4Step.hh"
61#include "G4Pow.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65G4double G4ICRU49NuclearStoppingModel::ad[] = {0.0};
66G4double G4ICRU49NuclearStoppingModel::ed[] = {0.0};
67
68using namespace std;
69
71 : G4VEmModel(nam),lossFlucFlag(false)
72{
73 theZieglerFactor = eV*cm2*1.0e-15;
74 g4pow = G4Pow::GetInstance();
75 if(ad[0] == 0.0) { InitialiseNuclearStopping(); }
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
81{}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86 const G4DataVector&)
87{}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
91void
93 std::vector<G4DynamicParticle*>*,
95 const G4DynamicParticle*,
97{}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
103 const G4Material* mat,
104 const G4ParticleDefinition* p,
105 G4double kinEnergy,
106 G4double)
107{
108 G4double nloss = 0.0;
109 if(kinEnergy <= 0.0) { return nloss; }
110
111 // projectile
112 G4double mass1 = p->GetPDGMass();
113 G4double z1 = std::fabs(p->GetPDGCharge()/eplus);
114
115 if(kinEnergy*proton_mass_c2/mass1 > z1*z1*MeV) { return nloss; }
116
117 // Projectile nucleus
118 mass1 /= amu_c2;
119
120 // loop for the elements in the material
121 G4int numberOfElements = mat->GetNumberOfElements();
122 const G4ElementVector* theElementVector = mat->GetElementVector();
123 const G4double* atomDensity = mat->GetAtomicNumDensityVector();
124
125 for (G4int iel=0; iel<numberOfElements; iel++) {
126 const G4Element* element = (*theElementVector)[iel] ;
127 G4double z2 = element->GetZ();
128 G4double mass2 = element->GetA()*mole/g ;
129 nloss += (NuclearStoppingPower(kinEnergy, z1, z2, mass1, mass2))
130 * atomDensity[iel] ;
131 }
132 nloss *= theZieglerFactor;
133 return nloss;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
139G4ICRU49NuclearStoppingModel::NuclearStoppingPower(G4double kineticEnergy,
140 G4double z1, G4double z2,
141 G4double mass1, G4double mass2)
142{
143 G4double energy = kineticEnergy/keV ; // energy in keV
144 G4double nloss = 0.0;
145 G4double z12 = z1*z2;
146 G4int iz1 = G4int(z1);
147 G4int iz2 = G4int(z2);
148
149 G4double rm;
150 if(iz1 > 1) { rm = (mass1 + mass2)*(g4pow->Z23(iz1) + g4pow->Z23(iz2)); }
151 else { rm = (mass1 + mass2)*g4pow->Z13(iz2); }
152
153 G4double er = 32.536 * mass2 * energy / ( z12 * rm ) ; // reduced energy
154
155 if (er >= ed[0]) { nloss = ad[0]; }
156 else {
157 // the table is inverse in energy
158 for (G4int i=102; i>=0; --i)
159 {
160 if (er <= ed[i]) {
161 nloss = (ad[i] - ad[i+1])*(er - ed[i+1])/(ed[i] - ed[i+1]) + ad[i+1];
162 break;
163 }
164 }
165 }
166
167 // Stragling
168 if(lossFlucFlag) {
169 // G4double sig = 4.0 * mass1 * mass2 / ((mass1 + mass2)*(mass1 + mass2)*
170 // (4.0 + 0.197*std::pow(er,-1.6991)+6.584*std::pow(er,-1.0494))) ;
171 G4double sig = 4.0 * mass1 * mass2 / ((mass1 + mass2)*(mass1 + mass2)*
172 (4.0 + 0.197/(er*er) + 6.584/er));
173
174 nloss *= G4RandGauss::shoot(1.0,sig);
175 lossFlucFlag = false;
176 }
177
178 nloss *= 8.462 * z12 * mass1 / rm; // Return to [ev/(10^15 atoms/cm^2]
179
180 if ( nloss < 0.0) { nloss = 0.0; }
181
182 return nloss;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
187void G4ICRU49NuclearStoppingModel::InitialiseNuclearStopping()
188{
189 const G4double nuca[104][2] = {
190 { 1.0E+8, 5.831E-8},
191 { 8.0E+7, 7.288E-8},
192 { 6.0E+7, 9.719E-8},
193 { 5.0E+7, 1.166E-7},
194 { 4.0E+7, 1.457E-7},
195 { 3.0E+7, 1.942E-7},
196 { 2.0E+7, 2.916E-7},
197 { 1.5E+7, 3.887E-7},
198
199 { 1.0E+7, 5.833E-7},
200 { 8.0E+6, 7.287E-7},
201 { 6.0E+6, 9.712E-7},
202 { 5.0E+6, 1.166E-6},
203 { 4.0E+6, 1.457E-6},
204 { 3.0E+6, 1.941E-6},
205 { 2.0E+6, 2.911E-6},
206 { 1.5E+6, 3.878E-6},
207
208 { 1.0E+6, 5.810E-6},
209 { 8.0E+5, 7.262E-6},
210 { 6.0E+5, 9.663E-6},
211 { 5.0E+5, 1.157E-5},
212 { 4.0E+5, 1.442E-5},
213 { 3.0E+5, 1.913E-5},
214 { 2.0E+5, 2.845E-5},
215 { 1.5E+5, 3.762E-5},
216
217 { 1.0E+5, 5.554E-5},
218 { 8.0E+4, 6.866E-5},
219 { 6.0E+4, 9.020E-5},
220 { 5.0E+4, 1.070E-4},
221 { 4.0E+4, 1.319E-4},
222 { 3.0E+4, 1.722E-4},
223 { 2.0E+4, 2.499E-4},
224 { 1.5E+4, 3.248E-4},
225
226 { 1.0E+4, 4.688E-4},
227 { 8.0E+3, 5.729E-4},
228 { 6.0E+3, 7.411E-4},
229 { 5.0E+3, 8.718E-4},
230 { 4.0E+3, 1.063E-3},
231 { 3.0E+3, 1.370E-3},
232 { 2.0E+3, 1.955E-3},
233 { 1.5E+3, 2.511E-3},
234
235 { 1.0E+3, 3.563E-3},
236 { 8.0E+2, 4.314E-3},
237 { 6.0E+2, 5.511E-3},
238 { 5.0E+2, 6.430E-3},
239 { 4.0E+2, 7.756E-3},
240 { 3.0E+2, 9.855E-3},
241 { 2.0E+2, 1.375E-2},
242 { 1.5E+2, 1.736E-2},
243
244 { 1.0E+2, 2.395E-2},
245 { 8.0E+1, 2.850E-2},
246 { 6.0E+1, 3.552E-2},
247 { 5.0E+1, 4.073E-2},
248 { 4.0E+1, 4.802E-2},
249 { 3.0E+1, 5.904E-2},
250 { 1.5E+1, 9.426E-2},
251
252 { 1.0E+1, 1.210E-1},
253 { 8.0E+0, 1.377E-1},
254 { 6.0E+0, 1.611E-1},
255 { 5.0E+0, 1.768E-1},
256 { 4.0E+0, 1.968E-1},
257 { 3.0E+0, 2.235E-1},
258 { 2.0E+0, 2.613E-1},
259 { 1.5E+0, 2.871E-1},
260
261 { 1.0E+0, 3.199E-1},
262 { 8.0E-1, 3.354E-1},
263 { 6.0E-1, 3.523E-1},
264 { 5.0E-1, 3.609E-1},
265 { 4.0E-1, 3.693E-1},
266 { 3.0E-1, 3.766E-1},
267 { 2.0E-1, 3.803E-1},
268 { 1.5E-1, 3.788E-1},
269
270 { 1.0E-1, 3.711E-1},
271 { 8.0E-2, 3.644E-1},
272 { 6.0E-2, 3.530E-1},
273 { 5.0E-2, 3.444E-1},
274 { 4.0E-2, 3.323E-1},
275 { 3.0E-2, 3.144E-1},
276 { 2.0E-2, 2.854E-1},
277 { 1.5E-2, 2.629E-1},
278
279 { 1.0E-2, 2.298E-1},
280 { 8.0E-3, 2.115E-1},
281 { 6.0E-3, 1.883E-1},
282 { 5.0E-3, 1.741E-1},
283 { 4.0E-3, 1.574E-1},
284 { 3.0E-3, 1.372E-1},
285 { 2.0E-3, 1.116E-1},
286 { 1.5E-3, 9.559E-2},
287
288 { 1.0E-3, 7.601E-2},
289 { 8.0E-4, 6.668E-2},
290 { 6.0E-4, 5.605E-2},
291 { 5.0E-4, 5.008E-2},
292 { 4.0E-4, 4.352E-2},
293 { 3.0E-4, 3.617E-2},
294 { 2.0E-4, 2.768E-2},
295 { 1.5E-4, 2.279E-2},
296
297 { 1.0E-4, 1.723E-2},
298 { 8.0E-5, 1.473E-2},
299 { 6.0E-5, 1.200E-2},
300 { 5.0E-5, 1.052E-2},
301 { 4.0E-5, 8.950E-3},
302 { 3.0E-5, 7.246E-3},
303 { 2.0E-5, 5.358E-3},
304 { 1.5E-5, 4.313E-3},
305 { 0.0, 3.166E-3}
306 };
307
308 for(G4int i=0; i<104; ++i) {
309 ed[i] = nuca[i][0];
310 ad[i] = nuca[i][1];
311 }
312}
313
314//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetA() const
Definition: G4Element.hh:138
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double)
G4ICRU49NuclearStoppingModel(const G4String &nam="ICRU49NucStopping")
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double Z23(G4int Z)
Definition: G4Pow.hh:134
G4double Z13(G4int Z)
Definition: G4Pow.hh:110