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