source: trunk/source/processes/electromagnetic/standard/src/G4ICRU49NuclearStoppingModel.cc@ 1199

Last change on this file since 1199 was 1197, checked in by garnier, 16 years ago

nvx fichiers dans CVS

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