Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ICRU73QOModel.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: G4ICRU73QOModel.cc 74581 2013-10-15 12:03:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4ICRU73QOModel
34 //
35 // Author: Alexander Bagulya
36 //
37 // Creation date: 21.05.2010
38 //
39 // Modifications:
40 //
41 //
42 // -------------------------------------------------------------------
43 //
44 
45 
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
49 #include "G4ICRU73QOModel.hh"
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "Randomize.hh"
53 #include "G4Electron.hh"
55 #include "G4LossTableManager.hh"
56 #include "G4AntiProton.hh"
57 #include "G4Log.hh"
58 #include "G4Exp.hh"
59 
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61 
62 using namespace std;
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67  : G4VEmModel(nam),
68  particle(0),
69  isInitialised(false)
70 {
71  mass = charge = chargeSquare = massRate = ratio = 0.0;
72  if(p) { SetParticle(p); }
73  SetHighEnergyLimit(10.0*MeV);
74 
75  lowestKinEnergy = 5.0*keV;
76 
77  sizeL0 = 67;
78  sizeL1 = 22;
79  sizeL2 = 14;
80 
81  theElectron = G4Electron::Electron();
82 
83  for (G4int i = 0; i < 100; ++i)
84  {
85  indexZ[i] = -1;
86  }
87  for(G4int i = 0; i < NQOELEM; ++i)
88  {
89  if(ZElementAvailable[i] > 0) {
90  indexZ[ZElementAvailable[i]] = i;
91  }
92  }
93  fParticleChange = 0;
94  denEffData = 0;
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98 
100 {}
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 
105  const G4DataVector&)
106 {
107  if(p != particle) SetParticle(p);
108 
109  // always false before the run
110  SetDeexcitationFlag(false);
111 
112  if(!isInitialised) {
113  isInitialised = true;
114 
115  G4String pname = particle->GetParticleName();
116  fParticleChange = GetParticleChangeForLoss();
118  denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
119  }
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 
125  const G4ParticleDefinition* p,
126  G4double kineticEnergy,
127  G4double cutEnergy,
128  G4double maxKinEnergy)
129 {
130  G4double cross = 0.0;
131  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
132  G4double maxEnergy = std::min(tmax,maxKinEnergy);
133  if(cutEnergy < maxEnergy) {
134 
135  G4double energy = kineticEnergy + mass;
136  G4double energy2 = energy*energy;
137  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
138  cross = 1.0/cutEnergy - 1.0/maxEnergy
139  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
140 
141  cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
142  }
143 
144  return cross;
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150  const G4ParticleDefinition* p,
151  G4double kineticEnergy,
152  G4double Z, G4double,
153  G4double cutEnergy,
154  G4double maxEnergy)
155 {
157  (p,kineticEnergy,cutEnergy,maxEnergy);
158  return cross;
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
164  const G4Material* material,
165  const G4ParticleDefinition* p,
166  G4double kineticEnergy,
167  G4double cutEnergy,
168  G4double maxEnergy)
169 {
170  G4double eDensity = material->GetElectronDensity();
172  (p,kineticEnergy,cutEnergy,maxEnergy);
173  return cross;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
179  const G4ParticleDefinition* p,
180  G4double kineticEnergy,
181  G4double cutEnergy)
182 {
183  SetParticle(p);
184  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
185  G4double tkin = kineticEnergy/massRate;
186  G4double dedx = 0.0;
187  if(tkin > lowestKinEnergy) { dedx = DEDX(material, tkin); }
188  else { dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy); }
189 
190  if (cutEnergy < tmax) {
191 
192  G4double tau = kineticEnergy/mass;
193  G4double gam = tau + 1.0;
194  G4double bg2 = tau * (tau+2.0);
195  G4double beta2 = bg2/(gam*gam);
196  G4double x = cutEnergy/tmax;
197 
198  dedx += chargeSquare*( G4Log(x) + (1.0 - x)*beta2 ) * twopi_mc2_rcl2
199  * material->GetElectronDensity()/beta2;
200  }
201  if(dedx < 0.0) { dedx = 0.0; }
202  return dedx;
203 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206 
207 G4double G4ICRU73QOModel::DEDX(const G4Material* material,
208  G4double kineticEnergy)
209 {
210  G4double eloss = 0.0;
211  const G4int numberOfElements = material->GetNumberOfElements();
212  const G4double* theAtomicNumDensityVector =
213  material->GetAtomicNumDensityVector();
214 
215  // Bragg's rule calculation
216  const G4ElementVector* theElementVector =
217  material->GetElementVector() ;
218 
219  // loop for the elements in the material
220  for (G4int i=0; i<numberOfElements; ++i)
221  {
222  const G4Element* element = (*theElementVector)[i] ;
223  eloss += DEDXPerElement(G4int(element->GetZ()), kineticEnergy)
224  * theAtomicNumDensityVector[i] * G4int(element->GetZ());
225  }
226  return eloss;
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230 
231 G4double G4ICRU73QOModel::DEDXPerElement(G4int AtomicNumber,
232  G4double kineticEnergy)
233 {
234  G4int Z = AtomicNumber;
235  if(Z > 97) { Z = 97; }
236  G4int nbOfShells = GetNumberOfShells(Z);
237  if(nbOfShells < 1) { nbOfShells = 1; }
238 
239  G4double v = CLHEP::c_light * std::sqrt( 2.0*kineticEnergy/proton_mass_c2 );
240 
241  G4double fBetheVelocity = CLHEP::fine_structure_const*CLHEP::c_light/v;
242 
243  G4double tau = kineticEnergy/proton_mass_c2;
244  G4double gam = tau + 1.0;
245  G4double bg2 = tau * (tau+2.0);
246  G4double beta2 = bg2/(gam*gam);
247 
248  G4double l0Term = 0, l1Term = 0, l2Term = 0;
249 
250  for (G4int nos = 0; nos < nbOfShells; ++nos){
251 
252  G4double NormalizedEnergy = (2.0*CLHEP::electron_mass_c2*beta2) /
253  GetShellEnergy(Z,nos);
254 
255  G4double shStrength = GetShellStrength(Z,nos);
256 
257  G4double l0 = GetL0(NormalizedEnergy);
258  l0Term += shStrength * l0;
259 
260  G4double l1 = GetL1(NormalizedEnergy);
261  l1Term += shStrength * l1;
262 
263  G4double l2 = GetL2(NormalizedEnergy);
264  l2Term += shStrength * l2;
265 
266  }
267  G4double dedx = 2*CLHEP::twopi_mc2_rcl2*chargeSquare*factorBethe[Z]*
268  (l0Term + charge*fBetheVelocity*l1Term
269  + chargeSquare*fBetheVelocity*fBetheVelocity*l2Term)/beta2;
270  return dedx;
271 }
272 
273 
274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
275 
276 G4double G4ICRU73QOModel::GetOscillatorEnergy(G4int Z,
277  G4int nbOfTheShell) const
278 {
279  G4int idx = denEffData->GetElementIndex(Z, kStateUndefined);
280  if(idx == -1) { idx = denEffData->GetElementIndex(Z-1, kStateUndefined); }
281  G4double PlasmaEnergy = denEffData->GetPlasmaEnergy(idx);
282 
283  G4double PlasmaEnergy2 = PlasmaEnergy * PlasmaEnergy;
284 
285  G4double plasmonTerm = 0.66667
286  * G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell)
287  * PlasmaEnergy2 / (Z*Z) ;
288 
289  G4double ionTerm = G4Exp(0.5) *
290  (G4AtomicShells::GetBindingEnergy(Z,nbOfTheShell)) ;
291  G4double ionTerm2 = ionTerm*ionTerm ;
292 
293  G4double oscShellEnergy = std::sqrt( ionTerm2 + plasmonTerm );
294 
295  return oscShellEnergy;
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
299 
300 G4double G4ICRU73QOModel::GetL0(G4double normEnergy) const
301 {
302  G4int n;
303 
304  for(n = 0; n < sizeL0; n++) {
305  if( normEnergy < L0[n][0] ) break;
306  }
307  if(0 == n) { n = 1; }
308  if(n >= sizeL0) { n = sizeL0 - 1; }
309 
310  G4double l0 = L0[n][1];
311  G4double l0p = L0[n-1][1];
312  G4double bethe = l0p + (l0 - l0p) * ( normEnergy - L0[n-1][0]) /
313  (L0[n][0] - L0[n-1][0]);
314 
315  return bethe ;
316 }
317 
318 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
319 
320 G4double G4ICRU73QOModel::GetL1(G4double normEnergy) const
321 {
322  G4int n;
323 
324  for(n = 0; n < sizeL1; n++) {
325  if( normEnergy < L1[n][0] ) break;
326  }
327  if(0 == n) n = 1 ;
328  if(n >= sizeL1) n = sizeL1 - 1 ;
329 
330  G4double l1 = L1[n][1];
331  G4double l1p = L1[n-1][1];
332  G4double barkas= l1p + (l1 - l1p) * ( normEnergy - L1[n-1][0]) /
333  (L1[n][0] - L1[n-1][0]);
334 
335  return barkas;
336 }
337 
338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
339 
340 G4double G4ICRU73QOModel::GetL2(G4double normEnergy) const
341 {
342  G4int n;
343  for(n = 0; n < sizeL2; n++) {
344  if( normEnergy < L2[n][0] ) break;
345  }
346  if(0 == n) n = 1 ;
347  if(n >= sizeL2) n = sizeL2 - 1 ;
348 
349  G4double l2 = L2[n][1];
350  G4double l2p = L2[n-1][1];
351  G4double bloch = l2p + (l2 - l2p) * ( normEnergy - L2[n-1][0]) /
352  (L2[n][0] - L2[n-1][0]);
353 
354  return bloch;
355 }
356 
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
358 
360  const G4DynamicParticle*,
361  G4double&,
362  G4double&,
363  G4double)
364 {}
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
367 
368 void G4ICRU73QOModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
369  const G4MaterialCutsCouple*,
370  const G4DynamicParticle* dp,
371  G4double xmin,
372  G4double maxEnergy)
373 {
374  G4double tmax = MaxSecondaryKinEnergy(dp);
375  G4double xmax = std::min(tmax, maxEnergy);
376  if(xmin >= xmax) { return; }
377 
378  G4double kineticEnergy = dp->GetKineticEnergy();
379  G4double energy = kineticEnergy + mass;
380  G4double energy2 = energy*energy;
381  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
382  G4double grej = 1.0;
383  G4double deltaKinEnergy, f;
384 
385  G4ThreeVector direction = dp->GetMomentumDirection();
386 
387  // sampling follows ...
388  do {
390  deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
391 
392  f = 1.0 - beta2*deltaKinEnergy/tmax;
393 
394  if(f > grej) {
395  G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
396  << "Majorant " << grej << " < "
397  << f << " for e= " << deltaKinEnergy
398  << G4endl;
399  }
400 
401  } while( grej*G4UniformRand() >= f );
402 
403  G4double deltaMomentum =
404  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
405  G4double totMomentum = energy*sqrt(beta2);
406  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
407  (deltaMomentum * totMomentum);
408  if(cost > 1.0) { cost = 1.0; }
409  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
410 
411  G4double phi = twopi * G4UniformRand() ;
412 
413  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
414  deltaDirection.rotateUz(direction);
415 
416  // Change kinematics of primary particle
417  kineticEnergy -= deltaKinEnergy;
418  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
419  finalP = finalP.unit();
420 
421  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
422  fParticleChange->SetProposedMomentumDirection(finalP);
423 
424  // create G4DynamicParticle object for delta ray
425  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,
426  deltaKinEnergy);
427 
428  vdp->push_back(delta);
429 }
430 
431 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
432 
434  G4double kinEnergy)
435 {
436  if(pd != particle) { SetParticle(pd); }
437  G4double tau = kinEnergy/mass;
438  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
439  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
440  return tmax;
441 }
442 
443 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
444 
445 const G4int G4ICRU73QOModel::ZElementAvailable[NQOELEM] = {1,2,4,6,7,8,10,13,14,-18,
446  22,26,28,29,32,36,42,47,
447  50,54,73,74,78,79,82,92};
448 
449 const G4int G4ICRU73QOModel::nbofShellsForElement[NQOELEM] = {1,1,2,3,3,3,3,4,5,4,
450  5,5,5,5,6,4,6,6,
451  7,6,6,8,7,7,9,9};
452 
453 const G4int G4ICRU73QOModel::startElemIndex[NQOELEM] = {0,1,2,4,7,10,13,16,20,25,
454  29,34,39,44,49,55,59,65,
455  71,78,84,90,98,105,112,121};
456 
457 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
458 
459 // SubShellOccupation = Z * ShellStrength
460 const G4double G4ICRU73QOModel::SubShellOccupation[NQODATA] =
461  {
462  1.000, // H 0
463  2.000, // He 1
464  1.930, 2.070, // Be 2-3
465  1.992, 1.841, 2.167, // C 4-6
466  1.741, 1.680, 3.579, // N 7-9
467  1.802, 1.849, 4.349, // O 10-12
468  1.788, 2.028, 6.184, // Ne 13-15
469  1.623, 2.147, 6.259, 2.971, // Al 16-19
470  1.631, 2.094, 6.588, 2.041, 1.646, // Si 20-24
471  1.535, 8.655, 1.706, 6.104, // Ar 25-28
472  1.581, 8.358, 8.183, 2.000, 1.878, // Ti 29-33
473  1.516, 8.325, 8.461, 6.579, 1.119, // Fe 34-38
474  1.422, 7.81, 8.385, 8.216, 2.167, // Ni 39-43
475  1.458, 8.049, 8.79, 9.695, 1.008, // Cu 44-48
476  1.442, 7.791, 7.837, 10.122, 2.463, 2.345, // Ge 49-54
477  1.645, 7.765, 19.192, 7.398, // Kr 55-58
478  1.313, 6.409, 19.229, 8.633, 5.036, 1.380, // Mo 59-64
479  1.295, 6.219, 18.751, 8.748, 10.184, 1.803, // Ag 65-70
480  1.277, 6.099, 20.386, 8.011, 10.007, 2.272, 1.948, // Sn 71-77
481  1.563, 6.312, 21.868, 5.762, 11.245, 7.250, // Xe 78-83
482  0.9198, 6.5408, 18.9727, 24.9149, 15.0161, 6.6284, // Ta 84-89
483  1.202, 5.582, 19.527, 18.741, 8.411, 14.387, 4.042, 2.108, // W 90-97
484  1.159, 5.467, 18.802, 33.905, 8.300, 9.342, 1.025, // Pt 98-104
485  1.124, 5.331, 18.078, 34.604, 8.127, 10.414, 1.322, // Au 105-111
486  2.000, 8.000, 18.000, 18.000, 14.000, 8.000, 10.000, 2.000, 2.000, // Pb 112-120
487  2.000, 8.000, 18.000, 32.000, 18.000, 8.000, 2.000, 1.000, 3.000 // U 121-129
488 };
489 
490 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
491 
492 // ShellEnergy in eV
493 const G4double G4ICRU73QOModel::ShellEnergy[NQODATA] =
494  {
495  19.2, // H
496  41.8, // He
497  209.11, 21.68, // Be
498  486.2, 60.95, 23.43, // C
499  732.61, 100.646, 23.550, // N
500  965.1, 129.85, 31.60, // O
501  1525.9, 234.9, 56.18, // Ne
502  2701, 476.5, 150.42, 16.89, // Al
503  3206.1, 586.4, 186.8, 23.52, 14.91, // Si
504  5551.6, 472.43, 124.85, 22.332, // Ar
505  8554.6, 850.58, 93.47, 39.19, 19.46, // Ti
506  12254.7, 1279.29, 200.35, 49.19, 17.66, // Fe
507  14346.9, 1532.28, 262.71, 74.37, 23.03, // Ni
508  15438.5, 1667.96, 294.1, 70.69, 16.447, // Cu
509  19022.1, 2150.79, 455.79, 179.87, 57.89, 20.95, // Ge
510  24643, 2906.4, 366.85, 22.24, // Kr
511  34394, 4365.3, 589.36, 129.42, 35.59, 18.42, // Mo
512  43664.3, 5824.91, 909.79, 175.47, 54.89, 19.63, // Ag
513  49948, 6818.2, 1036.1, 172.65, 70.89, 33.87, 14.54, // Sn
514  58987, 8159, 1296.6, 356.75, 101.03, 16.52, // Xe
515  88926, 18012, 3210, 575, 108.7, 30.8, // Ta
516  115025.9, 17827.44, 3214.36, 750.41, 305.21, 105.50, 38.09, 21.25, // W
517  128342, 20254, 3601.8, 608.1, 115.0, 42.75, 17.04, // Pt
518  131872, 20903, 3757.4, 682.1, 105.2, 44.89, 17.575, // Au
519  154449, 25067, 5105.0, 987.44, 247.59, 188.1, 40.61, 19.2, 15.17, // Pb
520  167282, 27868, 6022.7, 1020.4, 244.81, 51.33, 13, 11.06, 14.43 // U
521 };
522 
523 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
524 
525 // Data for L0 from: Sigmund P., Haagerup U. Phys. Rev. A34 (1986) 892-910
526 const G4double G4ICRU73QOModel::L0[67][2] =
527 {
528  {0.00, 0.000001},
529  {0.10, 0.000001},
530  {0.12, 0.00001},
531  {0.14, 0.00005},
532  {0.16, 0.00014},
533  {0.18, 0.00030},
534  {0.20, 0.00057},
535  {0.25, 0.00189},
536  {0.30, 0.00429},
537  {0.35, 0.00784},
538  {0.40, 0.01248},
539  {0.45, 0.01811},
540  {0.50, 0.02462},
541  {0.60, 0.03980},
542  {0.70, 0.05731},
543  {0.80, 0.07662},
544  {0.90, 0.09733},
545  {1.00, 0.11916},
546  {1.20, 0.16532},
547  {1.40, 0.21376},
548  {1.60, 0.26362},
549  {1.80, 0.31428},
550  {2.00, 0.36532},
551  {2.50, 0.49272},
552  {3.00, 0.61765},
553  {3.50, 0.73863},
554  {4.00, 0.85496},
555  {4.50, 0.96634},
556  {5.00, 1.07272},
557  {6.00, 1.27086},
558  {7.00, 1.45075},
559  {8.00, 1.61412},
560  {9.00, 1.76277},
561  {10.00, 1.89836},
562  {12.00, 2.13625},
563  {14.00, 2.33787},
564  {16.00, 2.51093},
565  {18.00, 2.66134},
566  {20.00, 2.79358},
567  {25.00, 3.06539},
568  {30.00, 3.27902},
569  {35.00, 3.45430},
570  {40.00, 3.60281},
571  {45.00, 3.73167},
572  {50.00, 3.84555},
573  {60.00, 4.04011},
574  {70.00, 4.20264},
575  {80.00, 4.34229},
576  {90.00, 4.46474},
577  {100.00, 4.57378},
578  {120.00, 4.76155},
579  {140.00, 4.91953},
580  {160.00, 5.05590},
581  {180.00, 5.17588},
582  {200.00, 5.28299},
583  {250.00, 5.50925},
584  {300.00, 5.69364},
585  {350.00, 5.84926},
586  {400.00, 5.98388},
587  {450.00, 6.10252},
588  {500.00, 6.20856},
589  {600.00, 6.39189},
590  {700.00, 6.54677},
591  {800.00, 6.68084},
592  {900.00, 6.79905},
593  {1000.00, 6.90474}
594 };
595 
596 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
597 
598 // Data for L1 from: Mikkelsen H.H., Sigmund P. Phys. Rev. A40 (1989) 101-116
599 const G4double G4ICRU73QOModel::L1[22][2] =
600 {
601  {0.00, -0.000001},
602  {0.10, -0.00001},
603  {0.20, -0.00049},
604  {0.30, -0.00084},
605  {0.40, 0.00085},
606  {0.50, 0.00519},
607  {0.60, 0.01198},
608  {0.70, 0.02074},
609  {0.80, 0.03133},
610  {0.90, 0.04369},
611  {1.00, 0.06035},
612  {2.00, 0.24023},
613  {3.00, 0.44284},
614  {4.00, 0.62012},
615  {5.00, 0.77031},
616  {6.00, 0.90390},
617  {7.00, 1.02705},
618  {8.00, 1.10867},
619  {9.00, 1.17546},
620  {10.00, 1.21599},
621  {15.00, 1.24349},
622  {20.00, 1.16752}
623 };
624 
625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
626 
627 // Data for L2 from: Mikkelsen H.H. Nucl. Instr. Meth. B58 (1991) 136-148
628 const G4double G4ICRU73QOModel::L2[14][2] =
629 {
630  {0.00, 0.000001},
631  {0.10, 0.00001},
632  {0.20, 0.00000},
633  {0.40, -0.00120},
634  {0.60, -0.00036},
635  {0.80, 0.00372},
636  {1.00, 0.01298},
637  {2.00, 0.08296},
638  {4.00, 0.21953},
639  {6.00, 0.23903},
640  {8.00, 0.20893},
641  {10.00, 0.10879},
642  {20.00, -0.88409},
643  {40.00, -1.13902}
644 };
645 
646 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
647 
648 // Correction obtained by V.Ivanchenko using G4BetheBlochModel
649 const G4double G4ICRU73QOModel::factorBethe[99] = { 1.0,
650 0.9637, 0.9872, 0.9469, 0.9875, 0.91, 0.989, 0.9507, 0.9773, 0.8621, 0.979, // 1 - 10
651 0.8357, 0.868, 0.9417, 0.9466, 0.8911, 0.905, 0.944, 0.9607, 0.928, 0.96, // 11 - 20
652 0.9098, 0.976, 0.8425, 0.8099, 0.7858, 0.947, 0.7248, 0.9106, 0.9246, 0.6821, // 21 - 30
653 0.7223, 0.9784, 0.774, 0.7953, 0.829, 0.9405, 0.8318, 0.8583, 0.8563, 0.8481, // 31 - 40
654 0.8207, 0.9033, 0.8063, 0.7837, 0.7818, 0.744, 0.875, 0.7693, 0.7871, 0.8459, // 41 - 50
655 0.8231, 0.8462, 0.853, 0.8736, 0.856, 0.8762, 0.8629, 0.8323, 0.8064, 0.7828, // 51 - 60
656 0.7533, 0.7273, 0.7093, 0.7157, 0.6823, 0.6612, 0.6418, 0.6395, 0.6323, 0.6221, // 61 - 70
657 0.6497, 0.6746, 0.8568, 0.8541, 0.6958, 0.6962, 0.7051, 0.863, 0.8588, 0.7226, // 71 - 80
658 0.7454, 0.78, 0.7783, 0.7996, 0.8216, 0.8632, 0.8558, 0.8792, 0.8745, 0.8676, // 81 - 90
659 0.8321, 0.8272, 0.7999, 0.7934, 0.7787, 0.7851, 0.7692, 0.7598};
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
G4int GetElementIndex(G4int Z, G4State mState)
const char * p
Definition: xmltok.h:285
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double GetZ() const
Definition: G4Element.hh:131
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
G4ICRU73QOModel(const G4ParticleDefinition *p=0, const G4String &nam="ICRU73QO")
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
string material
Definition: eplot.py:19
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:215
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
G4double GetPlasmaEnergy(G4int idx)
virtual ~G4ICRU73QOModel()
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
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
Hep3Vector unit() const
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double G4double
Definition: G4Types.hh:76
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)