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