source: trunk/source/processes/electromagnetic/standard/src/G4BraggIonModel.cc@ 1070

Last change on this file since 1070 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

File size: 24.2 KB
Line 
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: G4BraggIonModel.cc,v 1.24 2009/04/09 18:41:18 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4BraggIonModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 13.10.2004
39//
40// Modifications:
41// 11-05-05 Major optimisation of internal interfaces (V.Ivantchenko)
42// 29-11-05 Do not use G4Alpha class (V.Ivantchenko)
43// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
44// 25-04-06 Add stopping data from ASTAR (V.Ivanchenko)
45// 23-10-06 Reduce lowestKinEnergy to 0.25 keV (V.Ivanchenko)
46// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
47// CorrectionsAlongStep needed for ions(V.Ivanchenko)
48//
49
50// Class Description:
51//
52// Implementation of energy loss and delta-electron production by
53// slow charged heavy particles
54
55// -------------------------------------------------------------------
56//
57
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62#include "G4BraggIonModel.hh"
63#include "Randomize.hh"
64#include "G4Electron.hh"
65#include "G4ParticleChangeForLoss.hh"
66#include "G4LossTableManager.hh"
67#include "G4EmCorrections.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71using namespace std;
72
73G4BraggIonModel::G4BraggIonModel(const G4ParticleDefinition* p,
74 const G4String& nam)
75 : G4VEmModel(nam),
76 corr(0),
77 particle(0),
78 fParticleChange(0),
79 iMolecula(0),
80 isIon(false),
81 isInitialised(false)
82{
83 if(p) SetParticle(p);
84 SetHighEnergyLimit(2.0*MeV);
85
86 HeMass = 3.727417*GeV;
87 rateMassHe2p = HeMass/proton_mass_c2;
88 lowestKinEnergy = 1.0*keV/rateMassHe2p;
89 massFactor = 1000.*amu_c2/HeMass;
90 theZieglerFactor = eV*cm2*1.0e-15;
91 theElectron = G4Electron::Electron();
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
96G4BraggIonModel::~G4BraggIonModel()
97{}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101G4double G4BraggIonModel::MinEnergyCut(const G4ParticleDefinition*,
102 const G4MaterialCutsCouple* couple)
103{
104 return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109void G4BraggIonModel::Initialise(const G4ParticleDefinition* p,
110 const G4DataVector&)
111{
112 if(p != particle) SetParticle(p);
113
114 corrFactor = chargeSquare;
115
116 // always false before the run
117 SetDeexcitationFlag(false);
118
119 if(!isInitialised) {
120 isInitialised = true;
121
122 G4String pname = particle->GetParticleName();
123 if(particle->GetParticleType() == "nucleus" &&
124 pname != "deuteron" && pname != "triton") isIon = true;
125
126 corr = G4LossTableManager::Instance()->EmCorrections();
127
128 fParticleChange = GetParticleChangeForLoss();
129 }
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
134G4double G4BraggIonModel::GetChargeSquareRatio(const G4ParticleDefinition* p,
135 const G4Material* mat,
136 G4double kineticEnergy)
137{
138 // this method is called only for ions
139 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
140 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
141 return corrFactor;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146G4double G4BraggIonModel::GetParticleCharge(const G4ParticleDefinition* p,
147 const G4Material* mat,
148 G4double kineticEnergy)
149{
150 // this method is called only for ions
151 return corr->GetParticleCharge(p,mat,kineticEnergy);
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155
156G4double G4BraggIonModel::ComputeCrossSectionPerElectron(
157 const G4ParticleDefinition* p,
158 G4double kineticEnergy,
159 G4double cutEnergy,
160 G4double maxKinEnergy)
161{
162 G4double cross = 0.0;
163 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
164 G4double maxEnergy = std::min(tmax,maxKinEnergy);
165 if(cutEnergy < tmax) {
166
167 G4double energy = kineticEnergy + mass;
168 G4double energy2 = energy*energy;
169 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
170 cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax;
171
172 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
173 }
174 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
175 // << " tmax= " << tmax << " cross= " << cross << G4endl;
176
177 return cross;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
182G4double G4BraggIonModel::ComputeCrossSectionPerAtom(
183 const G4ParticleDefinition* p,
184 G4double kineticEnergy,
185 G4double Z, G4double,
186 G4double cutEnergy,
187 G4double maxEnergy)
188{
189 G4double cross = Z*ComputeCrossSectionPerElectron
190 (p,kineticEnergy,cutEnergy,maxEnergy);
191 return cross;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195
196G4double G4BraggIonModel::CrossSectionPerVolume(
197 const G4Material* material,
198 const G4ParticleDefinition* p,
199 G4double kineticEnergy,
200 G4double cutEnergy,
201 G4double maxEnergy)
202{
203 G4double eDensity = material->GetElectronDensity();
204 G4double cross = eDensity*ComputeCrossSectionPerElectron
205 (p,kineticEnergy,cutEnergy,maxEnergy);
206 return cross;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
211G4double G4BraggIonModel::ComputeDEDXPerVolume(const G4Material* material,
212 const G4ParticleDefinition* p,
213 G4double kineticEnergy,
214 G4double cutEnergy)
215{
216 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
217 G4double tmin = min(cutEnergy, tmax);
218 G4double tkin = kineticEnergy/massRate;
219 G4double dedx = 0.0;
220 if(tkin > lowestKinEnergy) dedx = DEDX(material, tkin);
221 else dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
222
223 if (cutEnergy < tmax) {
224
225 G4double tau = kineticEnergy/mass;
226 G4double gam = tau + 1.0;
227 G4double bg2 = tau * (tau+2.0);
228 G4double beta2 = bg2/(gam*gam);
229 G4double x = tmin/tmax;
230
231 dedx += (log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
232 * (material->GetElectronDensity())/beta2;
233 }
234
235 // now compute the total ionization loss
236
237 if (dedx < 0.0) dedx = 0.0 ;
238
239 dedx *= chargeSquare;
240
241 //G4cout << " tkin(MeV) = " << tkin/MeV << " dedx(MeVxcm^2/g) = "
242 // << dedx*gram/(MeV*cm2*material->GetDensity())
243 // << " q2 = " << chargeSquare << G4endl;
244
245 return dedx;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
249
250void G4BraggIonModel::CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
251 const G4DynamicParticle* dp,
252 G4double& eloss,
253 G4double&,
254 G4double length)
255{
256 // this method is called only for ions
257 const G4ParticleDefinition* p = dp->GetDefinition();
258 const G4Material* mat = couple->GetMaterial();
259 G4double preKinEnergy = dp->GetKineticEnergy();
260 G4double e = preKinEnergy - eloss*0.5;
261 if(e < 0.0) e = preKinEnergy*0.5;
262
263 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
264 GetModelOfFluctuations()->SetParticleAndCharge(p, q2);
265 eloss *= q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
266
267 if(nuclearStopping) {
268
269 G4double nloss = length*corr->NuclearDEDX(p,mat,e,false);
270
271 // too big energy loss
272 if(eloss + nloss > preKinEnergy) {
273 nloss *= (preKinEnergy/(eloss + nloss));
274 eloss = preKinEnergy;
275 } else {
276 eloss += nloss;
277 }
278 /*
279 G4cout << "G4ionIonisation::CorrectionsAlongStep: e= " << preKinEnergy
280 << " de= " << eloss << " NIEL= " << nloss
281 << " dynQ= " << dp->GetCharge()/eplus << G4endl;
282 */
283 fParticleChange->ProposeNonIonizingEnergyDeposit(nloss);
284 }
285}
286
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
288
289void G4BraggIonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
290 const G4MaterialCutsCouple*,
291 const G4DynamicParticle* dp,
292 G4double xmin,
293 G4double maxEnergy)
294{
295 G4double tmax = MaxSecondaryKinEnergy(dp);
296 G4double xmax = std::min(tmax, maxEnergy);
297 if(xmin >= xmax) return;
298
299 G4double kineticEnergy = dp->GetKineticEnergy();
300 G4double energy = kineticEnergy + mass;
301 G4double energy2 = energy*energy;
302 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
303 G4double grej = 1.0;
304 G4double deltaKinEnergy, f;
305
306 G4ThreeVector direction = dp->GetMomentumDirection();
307
308 // sampling follows ...
309 do {
310 G4double q = G4UniformRand();
311 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
312
313 f = 1.0 - beta2*deltaKinEnergy/tmax;
314
315 if(f > grej) {
316 G4cout << "G4BraggIonModel::SampleSecondary Warning! "
317 << "Majorant " << grej << " < "
318 << f << " for e= " << deltaKinEnergy
319 << G4endl;
320 }
321
322 } while( grej*G4UniformRand() >= f );
323
324 G4double deltaMomentum =
325 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
326 G4double totMomentum = energy*sqrt(beta2);
327 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
328 (deltaMomentum * totMomentum);
329 if(cost > 1.0) cost = 1.0;
330 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
331
332 G4double phi = twopi * G4UniformRand() ;
333
334 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
335 deltaDirection.rotateUz(direction);
336
337 // create G4DynamicParticle object for delta ray
338 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,
339 deltaKinEnergy);
340
341 vdp->push_back(delta);
342
343 // Change kinematics of primary particle
344 kineticEnergy -= deltaKinEnergy;
345 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
346 finalP = finalP.unit();
347
348 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
349 fParticleChange->SetProposedMomentumDirection(finalP);
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
353
354G4double G4BraggIonModel::MaxSecondaryEnergy(const G4ParticleDefinition* pd,
355 G4double kinEnergy)
356{
357 if(pd != particle) SetParticle(pd);
358 G4double tau = kinEnergy/mass;
359 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
360 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
361 return tmax;
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
365
366G4bool G4BraggIonModel::HasMaterial(const G4Material* material)
367{
368 const size_t numberOfMolecula = 11 ;
369 SetMoleculaNumber(numberOfMolecula) ;
370 G4String chFormula = material->GetChemicalFormula() ;
371
372 // ICRU Report N49, 1993. Ziegler model for He.
373 static G4String molName[numberOfMolecula] = {
374 "CaF_2", "Cellulose_Nitrate", "LiF", "Policarbonate",
375 "(C_2H_4)_N-Polyethylene", "(C_2H_4)_N-Polymethly_Methacralate",
376 "Polysterene", "SiO_2", "NaI", "H_2O",
377 "Graphite" } ;
378
379 // Search for the material in the table
380 for (size_t i=0; i<numberOfMolecula; i++) {
381 if (chFormula == molName[i]) {
382 SetMoleculaNumber(i) ;
383 return true ;
384 }
385 }
386 return false ;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
390
391G4double G4BraggIonModel::StoppingPower(const G4Material* material,
392 G4double kineticEnergy)
393{
394 G4double ionloss = 0.0 ;
395
396 if (iMolecula < 11) {
397
398 // The data and the fit from:
399 // ICRU Report N49, 1993. Ziegler's model for alpha
400 // He energy in internal units of parametrisation formula (MeV)
401
402 G4double T = kineticEnergy*rateMassHe2p/MeV ;
403
404 static G4double a[11][5] = {
405 {9.43672, 0.54398, 84.341, 1.3705, 57.422},
406 {67.1503, 0.41409, 404.512, 148.97, 20.99},
407 {5.11203, 0.453, 36.718, 50.6, 28.058},
408 {61.793, 0.48445, 361.537, 57.889, 50.674},
409 {7.83464, 0.49804, 160.452, 3.192, 0.71922},
410 {19.729, 0.52153, 162.341, 58.35, 25.668},
411 {26.4648, 0.50112, 188.913, 30.079, 16.509},
412 {7.8655, 0.5205, 63.96, 51.32, 67.775},
413 {8.8965, 0.5148, 339.36, 1.7205, 0.70423},
414 {2.959, 0.53255, 34.247, 60.655, 15.153},
415 {3.80133, 0.41590, 12.9966, 117.83, 242.28} };
416
417 static G4double atomicWeight[11] = {
418 101.96128, 44.0098, 16.0426, 28.0536, 42.0804,
419 104.1512, 44.665, 60.0843, 18.0152, 18.0152, 12.0};
420
421 G4int i = iMolecula;
422
423 // Free electron gas model
424 if ( T < 0.001 ) {
425 G4double slow = a[i][0] ;
426 G4double shigh = log( 1.0 + a[i][3]*1000.0 + a[i][4]*0.001 )
427 * a[i][2]*1000.0 ;
428 ionloss = slow*shigh / (slow + shigh) ;
429 ionloss *= sqrt(T*1000.0) ;
430
431 // Main parametrisation
432 } else {
433 G4double slow = a[i][0] * pow((T*1000.0), a[i][1]) ;
434 G4double shigh = log( 1.0 + a[i][3]/T + a[i][4]*T ) * a[i][2]/T ;
435 ionloss = slow*shigh / (slow + shigh) ;
436 /*
437 G4cout << "## " << i << ". T= " << T << " slow= " << slow
438 << " a0= " << a[i][0] << " a1= " << a[i][1]
439 << " shigh= " << shigh
440 << " dedx= " << ionloss << " q^2= " << HeEffChargeSquare(z, T*MeV)
441 << G4endl;
442 */
443 }
444 if ( ionloss < 0.0) ionloss = 0.0 ;
445
446 // He effective charge
447 G4double aa = atomicWeight[iMolecula];
448 ionloss /= (HeEffChargeSquare(0.5*aa, T)*aa);
449
450 // pure material (normally not the case for this function)
451 } else if(1 == (material->GetNumberOfElements())) {
452 G4double z = material->GetZ() ;
453 ionloss = ElectronicStoppingPower( z, kineticEnergy ) ;
454 }
455
456 return ionloss;
457}
458
459//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
460
461G4double G4BraggIonModel::ElectronicStoppingPower(G4double z,
462 G4double kineticEnergy) const
463{
464 G4double ionloss ;
465 G4int i = G4int(z)-1 ; // index of atom
466 if(i < 0) i = 0 ;
467 if(i > 91) i = 91 ;
468
469 // The data and the fit from:
470 // ICRU Report 49, 1993. Ziegler's type of parametrisations.
471 // Proton kinetic energy for parametrisation (keV/amu)
472
473 // He energy in internal units of parametrisation formula (MeV)
474 G4double T = kineticEnergy*rateMassHe2p/MeV ;
475
476 static G4double a[92][5] = {
477 {0.35485, 0.6456, 6.01525, 20.8933, 4.3515
478 },{ 0.58, 0.59, 6.3, 130.0, 44.07
479 },{ 1.42, 0.49, 12.25, 32.0, 9.161
480 },{ 2.206, 0.51, 15.32, 0.25, 8.995 //Be Ziegler77
481 // },{ 2.1895, 0.47183,7.2362, 134.30, 197.96 //Be from ICRU
482 },{ 3.691, 0.4128, 18.48, 50.72, 9.0
483 },{ 3.83523, 0.42993,12.6125, 227.41, 188.97
484 },{ 1.9259, 0.5550, 27.15125, 26.0665, 6.2768
485 },{ 2.81015, 0.4759, 50.0253, 10.556, 1.0382
486 },{ 1.533, 0.531, 40.44, 18.41, 2.718
487 },{ 2.303, 0.4861, 37.01, 37.96, 5.092
488 // Z= 11-20
489 },{ 9.894, 0.3081, 23.65, 0.384, 92.93
490 },{ 4.3, 0.47, 34.3, 3.3, 12.74
491 },{ 2.5, 0.625, 45.7, 0.1, 4.359
492 },{ 2.1, 0.65, 49.34, 1.788, 4.133
493 },{ 1.729, 0.6562, 53.41, 2.405, 3.845
494 },{ 1.402, 0.6791, 58.98, 3.528, 3.211
495 },{ 1.117, 0.7044, 69.69, 3.705, 2.156
496 },{ 2.291, 0.6284, 73.88, 4.478, 2.066
497 },{ 8.554, 0.3817, 83.61, 11.84, 1.875
498 },{ 6.297, 0.4622, 65.39, 10.14, 5.036
499 // Z= 21-30
500 },{ 5.307, 0.4918, 61.74, 12.4, 6.665
501 },{ 4.71, 0.5087, 65.28, 8.806, 5.948
502 },{ 6.151, 0.4524, 83.0, 18.31, 2.71
503 },{ 6.57, 0.4322, 84.76, 15.53, 2.779
504 },{ 5.738, 0.4492, 84.6, 14.18, 3.101
505 },{ 5.013, 0.4707, 85.8, 16.55, 3.211
506 },{ 4.32, 0.4947, 76.14, 10.85, 5.441
507 },{ 4.652, 0.4571, 80.73, 22.0, 4.952
508 },{ 3.114, 0.5236, 76.67, 7.62, 6.385
509 },{ 3.114, 0.5236, 76.67, 7.62, 7.502
510 // Z= 31-40
511 },{ 3.114, 0.5236, 76.67, 7.62, 8.514
512 },{ 5.746, 0.4662, 79.24, 1.185, 7.993
513 },{ 2.792, 0.6346, 106.1, 0.2986, 2.331
514 },{ 4.667, 0.5095, 124.3, 2.102, 1.667
515 },{ 2.44, 0.6346, 105.0, 0.83, 2.851
516 },{ 1.413, 0.7377, 147.9, 1.466, 1.016
517 },{ 11.72, 0.3826, 102.8, 9.231, 4.371
518 },{ 7.126, 0.4804, 119.3, 5.784, 2.454
519 },{ 11.61, 0.3955, 146.7, 7.031, 1.423
520 },{ 10.99, 0.41, 163.9, 7.1, 1.052
521 // Z= 41-50
522 },{ 9.241, 0.4275, 163.1, 7.954, 1.102
523 },{ 9.276, 0.418, 157.1, 8.038, 1.29
524 },{ 3.999, 0.6152, 97.6, 1.297, 5.792
525 },{ 4.306, 0.5658, 97.99, 5.514, 5.754
526 },{ 3.615, 0.6197, 86.26, 0.333, 8.689
527 },{ 5.8, 0.49, 147.2, 6.903, 1.289
528 },{ 5.6, 0.49, 130.0, 10.0, 2.844
529 },{ 3.55, 0.6068, 124.7, 1.112, 3.119
530 },{ 3.6, 0.62, 105.8, 0.1692, 6.026
531 },{ 5.4, 0.53, 103.1, 3.931, 7.767
532 // Z= 51-60
533 },{ 3.97, 0.6459, 131.8, 0.2233, 2.723
534 },{ 3.65, 0.64, 126.8, 0.6834, 3.411
535 },{ 3.118, 0.6519, 164.9, 1.208, 1.51
536 },{ 3.949, 0.6209, 200.5, 1.878, 0.9126
537 },{ 14.4, 0.3923, 152.5, 8.354, 2.597
538 },{ 10.99, 0.4599, 138.4, 4.811, 3.726
539 },{ 16.6, 0.3773, 224.1, 6.28, 0.9121
540 },{ 10.54, 0.4533, 159.3, 4.832, 2.529
541 },{ 10.33, 0.4502, 162.0, 5.132, 2.444
542 },{ 10.15, 0.4471, 165.6, 5.378, 2.328
543 // Z= 61-70
544 },{ 9.976, 0.4439, 168.0, 5.721, 2.258
545 },{ 9.804, 0.4408, 176.2, 5.675, 1.997
546 },{ 14.22, 0.363, 228.4, 7.024, 1.016
547 },{ 9.952, 0.4318, 233.5, 5.065, 0.9244
548 },{ 9.272, 0.4345, 210.0, 4.911, 1.258
549 },{ 10.13, 0.4146, 225.7, 5.525, 1.055
550 },{ 8.949, 0.4304, 213.3, 5.071, 1.221
551 },{ 11.94, 0.3783, 247.2, 6.655, 0.849
552 },{ 8.472, 0.4405, 195.5, 4.051, 1.604
553 },{ 8.301, 0.4399, 203.7, 3.667, 1.459
554 // Z= 71-80
555 },{ 6.567, 0.4858, 193.0, 2.65, 1.66
556 },{ 5.951, 0.5016, 196.1, 2.662, 1.589
557 },{ 7.495, 0.4523, 251.4, 3.433, 0.8619
558 },{ 6.335, 0.4825, 255.1, 2.834, 0.8228
559 },{ 4.314, 0.5558, 214.8, 2.354, 1.263
560 },{ 4.02, 0.5681, 219.9, 2.402, 1.191
561 },{ 3.836, 0.5765, 210.2, 2.742, 1.305
562 },{ 4.68, 0.5247, 244.7, 2.749, 0.8962
563 },{ 2.892, 0.6204, 208.6, 2.415, 1.416 //Au Z77
564 // },{ 3.223, 0.5883, 232.7, 2.954, 1.05 //Au ICRU
565 },{ 2.892, 0.6204, 208.6, 2.415, 1.416
566 // Z= 81-90
567 },{ 4.728, 0.5522, 217.0, 3.091, 1.386
568 },{ 6.18, 0.52, 170.0, 4.0, 3.224
569 },{ 9.0, 0.47, 198.0, 3.8, 2.032
570 },{ 2.324, 0.6997, 216.0, 1.599, 1.399
571 },{ 1.961, 0.7286, 223.0, 1.621, 1.296
572 },{ 1.75, 0.7427, 350.1, 0.9789, 0.5507
573 },{ 10.31, 0.4613, 261.2, 4.738, 0.9899
574 },{ 7.962, 0.519, 235.7, 4.347, 1.313
575 },{ 6.227, 0.5645, 231.9, 3.961, 1.379
576 },{ 5.246, 0.5947, 228.6, 4.027, 1.432
577 // Z= 91-92
578 },{ 5.408, 0.5811, 235.7, 3.961, 1.358
579 },{ 5.218, 0.5828, 245.0, 3.838, 1.25}
580 };
581
582 // Free electron gas model
583 if ( T < 0.001 ) {
584 G4double slow = a[i][0] ;
585 G4double shigh = log( 1.0 + a[i][3]*1000.0 + a[i][4]*0.001 )
586 * a[i][2]*1000.0 ;
587 ionloss = slow*shigh / (slow + shigh) ;
588 ionloss *= sqrt(T*1000.0) ;
589
590 // Main parametrisation
591 } else {
592 G4double slow = a[i][0] * pow((T*1000.0), a[i][1]) ;
593 G4double shigh = log( 1.0 + a[i][3]/T + a[i][4]*T ) * a[i][2]/T ;
594 ionloss = slow*shigh / (slow + shigh) ;
595 /*
596 G4cout << "## " << i << ". T= " << T << " slow= " << slow
597 << " a0= " << a[i][0] << " a1= " << a[i][1]
598 << " shigh= " << shigh
599 << " dedx= " << ionloss << " q^2= " << HeEffChargeSquare(z, T*MeV)
600 << G4endl;
601 */
602 }
603 if ( ionloss < 0.0) ionloss = 0.0 ;
604
605 // He effective charge
606 ionloss /= HeEffChargeSquare(z, T);
607
608 return ionloss;
609}
610
611//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
612
613G4double G4BraggIonModel::DEDX(const G4Material* material,
614 G4double kineticEnergy)
615{
616 G4double eloss = 0.0;
617 const G4int numberOfElements = material->GetNumberOfElements();
618 const G4double* theAtomicNumDensityVector =
619 material->GetAtomicNumDensityVector();
620
621 // compaund material with parametrisation
622 G4int iNist = astar.GetIndex(material);
623
624 if( iNist >= 0 ) {
625 G4double T = kineticEnergy*rateMassHe2p;
626 return astar.GetElectronicDEDX(iNist, T)*material->GetDensity()/
627 HeEffChargeSquare(astar.GetEffectiveZ(iNist), T/MeV);
628
629 } else if( HasMaterial(material) ) {
630
631 eloss = StoppingPower(material, kineticEnergy)*
632 material->GetDensity()/amu;
633
634 // pure material
635 } else if(1 == numberOfElements) {
636
637 G4double z = material->GetZ();
638 eloss = ElectronicStoppingPower(z, kineticEnergy)
639 * (material->GetTotNbOfAtomsPerVolume());
640
641 // Brugg's rule calculation
642 } else {
643 const G4ElementVector* theElementVector =
644 material->GetElementVector() ;
645
646 // loop for the elements in the material
647 for (G4int i=0; i<numberOfElements; i++)
648 {
649 const G4Element* element = (*theElementVector)[i] ;
650 eloss += ElectronicStoppingPower(element->GetZ(), kineticEnergy)
651 * theAtomicNumDensityVector[i];
652 }
653 }
654 return eloss*theZieglerFactor;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
658
659G4double G4BraggIonModel::HeEffChargeSquare(G4double z,
660 G4double kinEnergyHeInMeV) const
661{
662 // The aproximation of He effective charge from:
663 // J.F.Ziegler, J.P. Biersack, U. Littmark
664 // The Stopping and Range of Ions in Matter,
665 // Vol.1, Pergamon Press, 1985
666
667 static G4double c[6] = {0.2865, 0.1266, -0.001429,
668 0.02402,-0.01135, 0.001475};
669
670 G4double e = std::max(0.0,std::log(kinEnergyHeInMeV*massFactor));
671 G4double x = c[0] ;
672 G4double y = 1.0 ;
673 for (G4int i=1; i<6; i++) {
674 y *= e ;
675 x += y * c[i] ;
676 }
677
678 G4double w = 7.6 - e ;
679 w = 1.0 + (0.007 + 0.00005*z) * exp( -w*w ) ;
680 w = 4.0 * (1.0 - exp(-x)) * w * w ;
681
682 return w;
683}
684
685//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
686
Note: See TracBrowser for help on using the repository browser.