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