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

Last change on this file since 1330 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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