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

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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