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