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

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

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