Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmCorrections.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class
30//
31// File name: G4EmCorrections
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 13.01.2005
36//
37// Modifications:
38// 05.05.2005 V.Ivanchenko Fix misprint in Mott term
39// 26.11.2005 V.Ivanchenko Fix effective charge for heavy ions using original paper
40// 28.04.2006 V.Ivanchenko General cleanup, add finite size corrections
41// 13.05.2006 V.Ivanchenko Add corrections for ion stopping
42// 08.05.2007 V.Ivanchenko Use G4IonTable for ion mass instead of NistTable to avoid
43// division by zero
44// 29.02.2008 V.Ivanchenko use expantions for log and power function
45// 21.04.2008 Updated computations for ions (V.Ivanchenko)
46// 20.05.2008 Removed Finite Size correction (V.Ivanchenko)
47// 19.04.2012 Fix reproducibility problem (A.Ribon)
48//
49//
50// Class Description:
51//
52// This class provides calculation of EM corrections to ionisation
53//
54
55// -------------------------------------------------------------------
56//
57
58#include "G4EmCorrections.hh"
59#include "Randomize.hh"
61#include "G4SystemOfUnits.hh"
62#include "G4ParticleTable.hh"
63#include "G4IonTable.hh"
64#include "G4VEmModel.hh"
65#include "G4Proton.hh"
66#include "G4GenericIon.hh"
67#include "G4PhysicsLogVector.hh"
70#include "G4AtomicShells.hh"
71#include "G4Log.hh"
72#include "G4Exp.hh"
73#include "G4Pow.hh"
74#include "G4Threading.hh"
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
80
82 {0., 0., 0., 1.72, 2.09, 2.48, 2.82, 3.16, 3.53, 3.84, 4.15};
83const G4double G4EmCorrections::UK[20] = {1.9999, 2.0134, 2.0258, 2.0478, 2.0662,
84 2.0817, 2.0945, 2.0999, 2.1049, 2.1132,
85 2.1197, 2.1246, 2.1280, 2.1292, 2.1301,
86 2.1310, 2.1310, 2.1300, 2.1283, 2.1271};
87const G4double G4EmCorrections::VK[20] = {8.3410, 8.3373, 8.3340, 8.3287, 8.3247,
88 8.3219, 8.3201, 8.3194, 8.3191, 8.3188,
89 8.3191, 8.3199, 8.3211, 8.3218, 8.3226,
90 8.3244, 8.3264, 8.3285, 8.3308, 8.3320};
92const G4double G4EmCorrections::Eta[29] = {0.005, 0.007, 0.01, 0.015, 0.02,
93 0.03, 0.04, 0.05, 0.06, 0.08,
94 0.1, 0.15, 0.2, 0.3, 0.4,
95 0.5, 0.6, 0.7, 0.8, 1.0,
96 1.2, 1.4, 1.5, 1.7, 2.0, 3.0, 5.0, 7.0, 10.};
99const G4double G4EmCorrections::UL[] = {0.1215, 0.5265, 0.8411, 1.0878, 1.2828,
100 1.4379, 1.5032, 1.5617, 1.6608, 1.7401,
101 1.8036, 1.8543, 1.8756, 1.8945, 1.9262,
102 1.9508, 1.9696, 1.9836, 1.9890, 1.9935,
103 2.0001, 2.0039, 2.0053, 2.0049, 2.0040, 2.0028};
105
109
110#ifdef G4MULTITHREADED
111G4Mutex G4EmCorrections::theCorrMutex = G4MUTEX_INITIALIZER;
112#endif
113
115{
116 verbose = verb;
117 eth = 2.0*CLHEP::MeV;
118 eCorrMin = 25.*CLHEP::keV;
119 eCorrMax = 1.*CLHEP::GeV;
120
123
124 // G.S. Khandelwal Nucl. Phys. A116(1968)97 - 111.
125 // "Shell corrections for K- and L- electrons
126 nK = 20;
127 nL = 26;
128 nEtaK = 29;
129 nEtaL = 28;
130
131 // fill vectors
132 if(sBarkasCorr == nullptr) {
133#ifdef G4MULTITHREADED
134 G4MUTEXLOCK(&theCorrMutex);
135 if (sBarkasCorr == nullptr) {
136#endif
137 Initialise();
138 isMaster = true;
139#ifdef G4MULTITHREADED
140 }
141 G4MUTEXUNLOCK(&theCorrMutex);
142#endif
143 }
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
149{
150 for(G4int i=0; i<nIons; ++i) {delete stopData[i];}
151 if(isMaster) {
152 delete sBarkasCorr;
153 delete sThetaK;
154 delete sThetaL;
155 sBarkasCorr = sThetaK = sThetaL = nullptr;
156 }
157}
158
160 const G4Material* mat,
161 G4double kineticEnergy)
162{
163 if(kineticEnergy != kinEnergy || p != particle) {
164 particle = p;
165 kinEnergy = kineticEnergy;
166 mass = p->GetPDGMass();
167 tau = kineticEnergy / mass;
168 gamma = 1.0 + tau;
169 bg2 = tau * (tau+2.0);
170 beta2 = bg2/(gamma*gamma);
171 beta = std::sqrt(beta2);
172 ba2 = beta2/alpha2;
175 /(1. + 2.0*gamma*ratio + ratio*ratio);
177 if(charge > 1.5) { charge = effCharge.EffectiveCharge(p,mat,kinEnergy); }
178 q2 = charge*charge;
179 }
180 if(mat != material) {
181 material = mat;
185 }
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
191 const G4Material* mat,
193{
194 // . Z^3 Barkas effect in the stopping power of matter for charged particles
195 // J.C Ashley and R.H.Ritchie
196 // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
197 // and ICRU49 report
198 // valid for kineticEnergy < 0.5 MeV
199 // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
200
201 SetupKinematics(p, mat, e);
202 if(tau <= 0.0) { return 0.0; }
203
204 G4double Barkas = BarkasCorrection (p, mat, e);
205 G4double Bloch = BlochCorrection (p, mat, e);
206 G4double Mott = MottCorrection (p, mat, e);
207
208 G4double sum = (2.0*(Barkas + Bloch) + Mott);
209
210 if(verbose > 1) {
211 G4cout << "EmCorrections: E(MeV)= " << e/MeV << " Barkas= " << Barkas
212 << " Bloch= " << Bloch << " Mott= " << Mott
213 << " Sum= " << sum << " q2= " << q2 << G4endl;
214 G4cout << " ShellCorrection: " << ShellCorrection(p, mat, e)
215 << " Kshell= " << KShellCorrection(p, mat, e)
216 << " Lshell= " << LShellCorrection(p, mat, e)
217 << " " << mat->GetName() << G4endl;
218 }
220 return sum;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
226 const G4Material* mat,
227 G4double e)
228{
229 return 2.0*BarkasCorrection(p, mat, e)*
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234
236 const G4Material* mat,
237 G4double e)
238{
239 // . Z^3 Barkas effect in the stopping power of matter for charged particles
240 // J.C Ashley and R.H.Ritchie
241 // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
242 // and ICRU49 report
243 // valid for kineticEnergy < 0.5 MeV
244 // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
245 SetupKinematics(p, mat, e);
246 if(tau <= 0.0) { return 0.0; }
247
248 G4double Barkas = BarkasCorrection (p, mat, e);
249 G4double Bloch = BlochCorrection (p, mat, e);
250 G4double Mott = MottCorrection (p, mat, e);
251
252 G4double sum = 2.0*(Barkas*(charge - 1.0)/charge + Bloch) + Mott;
253
254 if(verbose > 1) {
255 G4cout << "EmCorrections: E(MeV)= " << e/MeV << " Barkas= " << Barkas
256 << " Bloch= " << Bloch << " Mott= " << Mott
257 << " Sum= " << sum << G4endl;
258 }
260
261 if(verbose > 1) { G4cout << " Sum= " << sum << G4endl; }
262 return sum;
263}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266
268 const G4MaterialCutsCouple* couple,
269 G4double e)
270{
271 // . Z^3 Barkas effect in the stopping power of matter for charged particles
272 // J.C Ashley and R.H.Ritchie
273 // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
274 // and ICRU49 report
275 // valid for kineticEnergy < 0.5 MeV
276 // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
277
278 G4double sum = 0.0;
279
280 if(ionHEModel) {
282 if(Z >= 100) Z = 99;
283 else if(Z < 1) Z = 1;
284
285 G4double ethscaled = eth*p->GetPDGMass()/proton_mass_c2;
286 G4int ionPDG = p->GetPDGEncoding();
287 if(thcorr.find(ionPDG)==thcorr.end()) { // Not found: fill the map
288 std::vector<G4double> v;
289 for(size_t i=0; i<ncouples; ++i){
290 v.push_back(ethscaled*ComputeIonCorrections(p,currmat[i],ethscaled));
291 }
292 thcorr.insert(std::pair< G4int, std::vector<G4double> >(ionPDG,v));
293 }
294
295 //G4cout << " map size=" << thcorr.size() << G4endl;
296 //for(std::map< G4int, std::vector<G4double> >::iterator
297 // it = thcorr.begin(); it != thcorr.end(); ++it){
298 // G4cout << "\t map element: first (key)=" << it->first
299 // << "\t second (vector): vec size=" << (it->second).size() << G4endl;
300 // for(size_t i=0; i<(it->second).size(); ++i){
301 // G4cout << "\t \t vec element: [" << i << "]=" << (it->second)[i]
302 //<< G4endl; } }
303
304 G4double rest = (thcorr.find(ionPDG)->second)[couple->GetIndex()];
305
306 sum = ComputeIonCorrections(p,couple->GetMaterial(),e) - rest/e;
307
308 if(verbose > 1) {
309 G4cout << " Sum= " << sum << " dSum= " << rest/e << G4endl;
310 }
311 }
312 return sum;
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316
318 const G4Material* mat,
319 G4double e)
320{
321 SetupKinematics(p, mat, e);
323 const G4double eexc2 = eexc*eexc;
324 return 0.5*G4Log(2.0*electron_mass_c2*bg2*tmax/eexc2)-beta2;
325}
326
327//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
328
330 const G4Material* mat,
331 G4double e)
332{
333 SetupKinematics(p, mat, e);
334 const G4double dedx = 0.5*tmax/(kinEnergy + mass);
335 return 0.5*dedx*dedx;
336}
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
339
341 const G4Material* mat,
342 G4double e)
343{
344 SetupKinematics(p, mat, e);
345 G4double term = 0.0;
346 for (G4int i = 0; i<numberOfElements; ++i) {
347
348 G4double Z = (*theElementVector)[i]->GetZ();
349 G4int iz = (*theElementVector)[i]->GetZasInt();
350 G4double f = 1.0;
351 G4double Z2= (Z-0.3)*(Z-0.3);
352 if(1 == iz) {
353 f = 0.5;
354 Z2 = 1.0;
355 }
356 G4double eta = ba2/Z2;
357 G4double tet = Z2*(1. + Z2*0.25*alpha2);
358 if(11 < iz) { tet = sThetaK->Value(Z); }
359 term += f*atomDensity[i]*KShell(tet,eta)/Z;
360 }
361
363
364 return term;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
368
370 const G4Material* mat,
371 G4double e)
372{
373 SetupKinematics(p, mat, e);
374 G4double term = 0.0;
375 for (G4int i = 0; i<numberOfElements; ++i) {
376
377 G4double Z = (*theElementVector)[i]->GetZ();
378 G4int iz = (*theElementVector)[i]->GetZasInt();
379 if(2 < iz) {
380 G4double Zeff = Z - ZD[10];
381 if(iz < 10) { Zeff = Z - ZD[iz]; }
382 G4double Z2= Zeff*Zeff;
383 G4double f = 0.125;
384 G4double eta = ba2/Z2;
385 G4double tet = sThetaL->Value(Z);
387 for(G4int j=1; j<nmax; ++j) {
389 if(15 >= iz) {
390 if(3 > j) { tet = 0.25*Z2*(1.0 + 5*Z2*alpha2/16.); }
391 else { tet = 0.25*Z2*(1.0 + Z2*alpha2/16.); }
392 }
393 //G4cout << " LShell: j= " << j << " ne= " << ne << " e(eV)= " << e/eV
394 // << " ThetaL= " << tet << G4endl;
395 term += f*ne*atomDensity[i]*LShell(tet,eta)/Z;
396 }
397 }
398 }
399
401
402 return term;
403}
404
405//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
406
408{
409 G4double corr = 0.0;
410
411 static const G4double TheK[20] =
412 {0.64, 0.65, 0.66, 0.68, 0.70, 0.72, 0.74, 0.75, 0.76, 0.78,
413 0.80, 0.82, 0.84, 0.85, 0.86, 0.88, 0.90, 0.92, 0.94, 0.95};
414
415
416 G4double x = tet;
417 G4int itet = 0;
418 G4int ieta = 0;
419 if(tet < TheK[0]) {
420 x = TheK[0];
421 } else if(tet > TheK[nK-1]) {
422 x = TheK[nK-1];
423 itet = nK-2;
424 } else {
425 itet = Index(x, TheK, nK);
426 }
427 // assimptotic case
428 if(eta >= Eta[nEtaK-1]) {
429 corr =
430 (Value(x, TheK[itet], TheK[itet+1], UK[itet], UK[itet+1]) +
431 Value(x, TheK[itet], TheK[itet+1], VK[itet], VK[itet+1])/eta +
432 Value(x, TheK[itet], TheK[itet+1], ZK[itet], ZK[itet+1])/(eta*eta))/eta;
433 } else {
434 G4double y = eta;
435 if(eta < Eta[0]) {
436 y = Eta[0];
437 } else {
438 ieta = Index(y, Eta, nEtaK);
439 }
440 corr = Value2(x, y, TheK[itet], TheK[itet+1], Eta[ieta], Eta[ieta+1],
441 CK[itet][ieta], CK[itet+1][ieta],
442 CK[itet][ieta+1], CK[itet+1][ieta+1]);
443 //G4cout << " x= " <<x<<" y= "<<y<<" tet= " <<TheK[itet]
444 // <<" "<< TheK[itet+1]<<" eta= "<< Eta[ieta]<<" "<< Eta[ieta+1]
445 // <<" CK= " << CK[itet][ieta]<<" "<< CK[itet+1][ieta]
446 // <<" "<< CK[itet][ieta+1]<<" "<< CK[itet+1][ieta+1]<<G4endl;
447 }
448 //G4cout << "Kshell: tet= " << tet << " eta= " << eta << " C= " << corr
449 // << " itet= " << itet << " ieta= " << ieta <<G4endl;
450 return corr;
451}
452
453//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
454
456{
457 G4double corr = 0.0;
458
459 static const G4double TheL[26] =
460 {0.24, 0.26, 0.28, 0.30, 0.32, 0.34, 0.35, 0.36, 0.38, 0.40,
461 0.42, 0.44, 0.45, 0.46, 0.48, 0.50, 0.52, 0.54, 0.55, 0.56,
462 0.58, 0.60, 0.62, 0.64, 0.65, 0.66};
463
464 G4double x = tet;
465 G4int itet = 0;
466 G4int ieta = 0;
467 if(tet < TheL[0]) {
468 x = TheL[0];
469 } else if(tet > TheL[nL-1]) {
470 x = TheL[nL-1];
471 itet = nL-2;
472 } else {
473 itet = Index(x, TheL, nL);
474 }
475
476 // assimptotic case
477 if(eta >= Eta[nEtaL-1]) {
478 corr = (Value(x, TheL[itet], TheL[itet+1], UL[itet], UL[itet+1])
479 + Value(x, TheL[itet], TheL[itet+1], VL[itet], VL[itet+1])/eta
480 )/eta;
481 } else {
482 G4double y = eta;
483 if(eta < Eta[0]) {
484 y = Eta[0];
485 } else {
486 ieta = Index(y, Eta, nEtaL);
487 }
488 corr = Value2(x, y, TheL[itet], TheL[itet+1], Eta[ieta], Eta[ieta+1],
489 CL[itet][ieta], CL[itet+1][ieta],
490 CL[itet][ieta+1], CL[itet+1][ieta+1]);
491 //G4cout << " x= " <<x<<" y= "<<y<<" tet= " <<TheL[itet]
492 // <<" "<< TheL[itet+1]<<" eta= "<< Eta[ieta]<<" "<< Eta[ieta+1]
493 // <<" CL= " << CL[itet][ieta]<<" "<< CL[itet+1][ieta]
494 // <<" "<< CL[itet][ieta+1]<<" "<< CL[itet+1][ieta+1]<<G4endl;
495 }
496 //G4cout<<"Lshell: tet= "<<tet<<" eta= "<<eta<<" itet= "<<itet
497 // <<" ieta= "<<ieta<<" Corr= "<<corr<<G4endl;
498 return corr;
499}
500
501//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
502
504 const G4Material* mat,
505 G4double e)
506{
507 SetupKinematics(p, mat, e);
508 G4double taulim= 8.0*MeV/mass;
509 G4double bg2lim= taulim * (taulim+2.0);
510
511 G4double* shellCorrectionVector =
513 G4double sh = 0.0;
514 G4double x = 1.0;
516
517 if ( bg2 >= bg2lim ) {
518 for (G4int k=0; k<3; ++k) {
519 x *= bg2 ;
520 sh += shellCorrectionVector[k]/x;
521 }
522
523 } else {
524 for (G4int k=0; k<3; ++k) {
525 x *= bg2lim ;
526 sh += shellCorrectionVector[k]/x;
527 }
528 sh *= G4Log(tau/taul)/G4Log(taulim/taul);
529 }
530 sh *= 0.5;
531 return sh;
532}
533
534
535//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
536
538 const G4Material* mat,
539 G4double ekin)
540{
541 SetupKinematics(p, mat, ekin);
542
543 G4double term = 0.0;
544 //G4cout << "### G4EmCorrections::ShellCorrection " << mat->GetName()
545 // << " " << ekin/MeV << " MeV " << G4endl;
546 for (G4int i = 0; i<numberOfElements; ++i) {
547
548 G4double res = 0.0;
549 G4double res0 = 0.0;
550 G4double Z = (*theElementVector)[i]->GetZ();
551 G4int iz = (*theElementVector)[i]->GetZasInt();
552 G4double Z2= (Z-0.3)*(Z-0.3);
553 G4double f = 1.0;
554 if(1 == iz) {
555 f = 0.5;
556 Z2 = 1.0;
557 }
558 G4double eta = ba2/Z2;
559 G4double tet = Z2*(1. + Z2*0.25*alpha2);
560 if(11 < iz) { tet = sThetaK->Value(Z); }
561 res0 = f*KShell(tet,eta);
562 res += res0;
563 //G4cout << " Z= " << iz << " Shell 0" << " tet= " << tet
564 // << " eta= " << eta << " resK= " << res0 << G4endl;
565 if(2 < iz) {
566 G4double Zeff = Z - ZD[10];
567 if(iz < 10) { Zeff = Z - ZD[iz]; }
568 Z2= Zeff*Zeff;
569 eta = ba2/Z2;
570 f = 0.125;
571 tet = sThetaL->Value(Z);
573 G4int nmax = std::min(4, ntot);
574 G4double norm = 0.0;
575 G4double eshell = 0.0;
576 for(G4int j=1; j<nmax; ++j) {
578 if(15 >= iz) {
579 if(3 > j) { tet = 0.25*Z2*(1.0 + 5*Z2*alpha2/16.); }
580 else { tet = 0.25*Z2*(1.0 + Z2*alpha2/16.); }
581 }
582 norm += ne;
583 eshell += tet*ne;
584 res0 = f*ne*LShell(tet,eta);
585 res += res0;
586 //G4cout << " Z= " << iz << " Shell " << j << " Ne= " << ne
587 // << " tet= " << tet << " eta= " << eta
588 // << " resL= " << res0 << G4endl;
589 }
590 if(ntot > nmax) {
591 eshell /= norm;
592
593 static const G4double HM[53] = {
594 12.0, 12.0, 12.0, 12.0, 11.9, 11.7, 11.5, 11.2, 10.8, 10.4,
595 10.0, 9.51, 8.97, 8.52, 8.03, 7.46, 6.95, 6.53, 6.18, 5.87,
596 5.61, 5.39, 5.19, 5.01, 4.86, 4.72, 4.62, 4.53, 4.44, 4.38,
597 4.32, 4.26, 4.20, 4.15, 4.1, 4.04, 4.00, 3.95, 3.93, 3.91,
598 3.90, 3.89, 3.89, 3.88, 3.88, 3.88, 3.88, 3.88, 3.89, 3.89,
599 3.90, 3.92, 3.93 };
600 static const G4double HN[31] = {
601 75.5, 61.9, 52.2, 45.1, 39.6, 35.4, 31.9, 29.1, 27.2, 25.8,
602 24.5, 23.6, 22.7, 22.0, 21.4, 20.9, 20.5, 20.2, 19.9, 19.7,
603 19.5, 19.3, 19.2, 19.1, 18.4, 18.8, 18.7, 18.6, 18.5, 18.4,
604 18.2};
605
606 // Add M-shell
607 if(28 > iz) {
608 res += f*(iz - 10)*LShell(eshell,HM[iz-11]*eta);
609 } else if(63 > iz) {
610 res += f*18*LShell(eshell,HM[iz-11]*eta);
611 } else {
612 res += f*18*LShell(eshell,HM[52]*eta);
613 }
614 // Add N-shell
615 if(32 < iz) {
616 if(60 > iz) {
617 res += f*(iz - 28)*LShell(eshell,HN[iz-33]*eta);
618 } else if(63 > iz) {
619 res += 4*LShell(eshell,HN[iz-33]*eta);
620 } else {
621 res += 4*LShell(eshell,HN[30]*eta);
622 }
623 // Add O-P-shells
624 if(60 < iz) {
625 res += f*(iz - 60)*LShell(eshell,150*eta);
626 }
627 }
628 }
629 }
630 term += res*atomDensity[i]/Z;
631 }
632
634 //G4cout << "# Shell Correction= " << term << G4endl;
635 return term;
636}
637
638//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
639
641 const G4Material* mat,
642 G4double e)
643{
644 SetupKinematics(p, mat, e);
645
651
652 G4double dedx = 0.0;
653
654 // density correction
655 static const G4double twoln10 = 2.0*G4Log(10.0);
656 G4double x = G4Log(bg2)/twoln10;
657 if ( x >= x0den ) {
658 dedx = twoln10*x - cden ;
659 if ( x < x1den ) dedx += aden*G4Exp(G4Log(x1den-x)*mden) ;
660 }
661
662 return dedx;
663}
664
665//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
666
668 const G4Material* mat,
669 G4double e)
670{
671 // . Z^3 Barkas effect in the stopping power of matter for charged particles
672 // J.C Ashley and R.H.Ritchie
673 // Physical review B Vol.5 No.7 1 April 1972 pp. 2393-2397
674 // valid for kineticEnergy > 0.5 MeV
675
676 SetupKinematics(p, mat, e);
677 G4double BarkasTerm = 0.0;
678
679 for (G4int i = 0; i<numberOfElements; ++i) {
680
681 G4double Z = (*theElementVector)[i]->GetZ();
682 G4int iz = (*theElementVector)[i]->GetZasInt();
683 if(iz == 47) {
684 BarkasTerm += atomDensity[i]*0.006812*G4Exp(-G4Log(beta)*0.9);
685 } else if(iz >= 64) {
686 BarkasTerm += atomDensity[i]*0.002833*G4Exp(-G4Log(beta)*1.2);
687 } else {
688
689 G4double X = ba2 / Z;
690 G4double b = 1.3;
691 if(1 == iz) {
692 if(material->GetName() == "G4_lH2") { b = 0.6; }
693 else { b = 1.8; }
694 }
695 else if(2 == iz) { b = 0.6; }
696 else if(10 >= iz) { b = 1.8; }
697 else if(17 >= iz) { b = 1.4; }
698 else if(18 == iz) { b = 1.8; }
699 else if(25 >= iz) { b = 1.4; }
700 else if(50 >= iz) { b = 1.35;}
701
702 G4double W = b/std::sqrt(X);
703
704 G4double val = sBarkasCorr->Value(W);
705 if(W > sBarkasCorr->Energy(46)) {
706 val *= sBarkasCorr->Energy(46)/W;
707 }
708 // G4cout << "i= " << i << " b= " << b << " W= " << W
709 // << " Z= " << Z << " X= " << X << " val= " << val<< G4endl;
710 BarkasTerm += val*atomDensity[i] / (std::sqrt(Z*X)*X);
711 }
712 }
713
714 BarkasTerm *= 1.29*charge/material->GetTotNbOfAtomsPerVolume();
715
716 return BarkasTerm;
717}
718
719//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
720
722 const G4Material* mat,
723 G4double e)
724{
725 SetupKinematics(p, mat, e);
726
727 G4double y2 = q2/ba2;
728
729 G4double term = 1.0/(1.0 + y2);
730 G4double del;
731 G4double j = 1.0;
732 do {
733 j += 1.0;
734 del = 1.0/(j* (j*j + y2));
735 term += del;
736 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
737 } while (del > 0.01*term);
738
739 G4double res = -y2*term;
740 return res;
741}
742
743//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
744
746 const G4Material* mat,
747 G4double e)
748{
749 SetupKinematics(p, mat, e);
751 return mterm;
752}
753
754//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
755
758 const G4Material* mat,
759 G4double ekin)
760{
761 G4double factor = 1.0;
762 if(p->GetPDGCharge() <= 2.5*CLHEP::eplus || nIons <= 0) { return factor; }
763
764 if(verbose > 1) {
765 G4cout << "EffectiveChargeCorrection: " << p->GetParticleName()
766 << " in " << mat->GetName()
767 << " ekin(MeV)= " << ekin/MeV << G4endl;
768 }
769
770 if(p != curParticle || mat != curMaterial) {
771 curParticle = p;
772 curMaterial = mat;
773 curVector = nullptr;
775 if(verbose > 1) {
776 G4cout << "G4EmCorrections::EffectiveChargeCorrection: Zion= "
777 << currentZ << " Aion= " << p->GetPDGMass()/amu_c2 << G4endl;
778 }
780 idx = -1;
781
782 for(G4int i=0; i<nIons; ++i) {
783 if(materialList[i] == mat && currentZ == Zion[i]) {
784 idx = i;
785 break;
786 }
787 }
788 //G4cout << " idx= " << idx << " dz= " << G4endl;
789 if(idx >= 0) {
790 if(nullptr == ionList[idx]) { BuildCorrectionVector(); }
792 } else {
793 return factor;
794 }
795 }
796 if(nullptr != curVector) {
797 factor = curVector->Value(ekin*massFactor);
798 if(verbose > 1) {
799 G4cout << "E= " << ekin << " factor= " << factor << " massfactor= "
800 << massFactor << G4endl;
801 }
802 }
803 return factor;
804}
805
806//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807
809 const G4String& mname,
810 G4PhysicsVector* dVector)
811{
812 G4int i = 0;
813 for(; i<nIons; ++i) {
814 if(Z == Zion[i] && A == Aion[i] && mname == materialName[i]) break;
815 }
816 if(i == nIons) {
817 Zion.push_back(Z);
818 Aion.push_back(A);
819 materialName.push_back(mname);
820 materialList.push_back(nullptr);
821 ionList.push_back(nullptr);
822 stopData.push_back(dVector);
823 nIons++;
824 if(verbose > 1) {
825 G4cout << "AddStoppingData Z= " << Z << " A= " << A << " " << mname
826 << " idx= " << i << G4endl;
827 }
828 }
829}
830
831//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
832
834{
835 if(nullptr == ionLEModel || nullptr == ionHEModel) {
836 return;
837 }
838
841 G4int Z = Zion[idx];
842 G4double A = Aion[idx];
844
845 if(verbose > 1) {
846 G4cout << "BuildCorrectionVector: Stopping for "
847 << curParticle->GetParticleName() << " in "
848 << materialName[idx] << " Ion Z= " << Z << " A= " << A
849 << " massFactor= " << massFactor << G4endl;
850 G4cout << " Nbins=" << nbinCorr << " Emin(MeV)=" << eCorrMin
851 << " Emax(MeV)=" << eCorrMax << " ion: "
852 << ion->GetParticleName() << G4endl;
853 }
854
855 G4PhysicsLogVector* vv =
857 const G4double eth0 = v->Energy(0);
858 const G4double escal = eth/massFactor;
859 G4double qe =
861 const G4double dedxt =
863 const G4double dedx1t =
866 const G4double rest = escal*(dedxt - dedx1t);
867 if(verbose > 1) {
868 G4cout << "Escal(MeV)= " << escal << " qe=" << qe
869 << " dedxt= " << dedxt << " dedx1t= " << dedx1t << G4endl;
870 }
871 for(G4int i=0; i<=nbinCorr; ++i) {
872 // energy in the local table (GenericIon)
873 G4double e = vv->Energy(i);
874 // energy of the real ion
875 G4double eion = e/massFactor;
876 // energy in the imput stopping data vector
877 G4double e1 = eion/A;
878 G4double dedx = (e1 < eth0)
879 ? v->Value(eth0)*std::sqrt(e1/eth0) : v->Value(e1);
881 G4double dedx1 = (e <= eth)
883 : ionHEModel->ComputeDEDXPerVolume(curMaterial, gion, e, e)*qe +
885 vv->PutValue(i, dedx/dedx1);
886 if(verbose > 1) {
887 G4cout << "E(MeV)=" << e/CLHEP::MeV << " Eion=" << eion/CLHEP::MeV
888 << " e1=" << e1 << " dedxRatio= " << dedx/dedx1
889 << " dedx=" << dedx << " dedx1=" << dedx1
890 << " qe=" << qe << " rest/eion=" << rest/eion << G4endl;
891 }
892 }
893 delete v;
894 ionList[idx] = ion;
895 stopData[idx] = vv;
896 if(verbose > 1) { G4cout << "End data set " << G4endl; }
897}
898
899//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
900
902{
904 ncouples = tb->GetTableSize();
905 if(currmat.size() != ncouples) {
906 currmat.resize(ncouples);
907 for(std::map< G4int, std::vector<G4double> >::iterator it =
908 thcorr.begin(); it != thcorr.end(); ++it){
909 (it->second).clear();
910 }
911 thcorr.clear();
912 for(size_t i=0; i<ncouples; ++i) {
914 G4String nam = currmat[i]->GetName();
915 for(G4int j=0; j<nIons; ++j) {
916 if(nam == materialName[j]) { materialList[j] = currmat[i]; }
917 }
918 }
919 }
920}
921
922//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
923
925{
926 // Z^3 Barkas effect in the stopping power of matter for charged particles
927 // J.C Ashley and R.H.Ritchie
928 // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
929 G4int i, j;
930 static const G4double fTable[47][2] = {
931 { 0.02, 21.5},
932 { 0.03, 20.0},
933 { 0.04, 18.0},
934 { 0.05, 15.6},
935 { 0.06, 15.0},
936 { 0.07, 14.0},
937 { 0.08, 13.5},
938 { 0.09, 13.},
939 { 0.1, 12.2},
940 { 0.2, 9.25},
941 { 0.3, 7.0},
942 { 0.4, 6.0},
943 { 0.5, 4.5},
944 { 0.6, 3.5},
945 { 0.7, 3.0},
946 { 0.8, 2.5},
947 { 0.9, 2.0},
948 { 1.0, 1.7},
949 { 1.2, 1.2},
950 { 1.3, 1.0},
951 { 1.4, 0.86},
952 { 1.5, 0.7},
953 { 1.6, 0.61},
954 { 1.7, 0.52},
955 { 1.8, 0.5},
956 { 1.9, 0.43},
957 { 2.0, 0.42},
958 { 2.1, 0.3},
959 { 2.4, 0.2},
960 { 3.0, 0.13},
961 { 3.08, 0.1},
962 { 3.1, 0.09},
963 { 3.3, 0.08},
964 { 3.5, 0.07},
965 { 3.8, 0.06},
966 { 4.0, 0.051},
967 { 4.1, 0.04},
968 { 4.8, 0.03},
969 { 5.0, 0.024},
970 { 5.1, 0.02},
971 { 6.0, 0.013},
972 { 6.5, 0.01},
973 { 7.0, 0.009},
974 { 7.1, 0.008},
975 { 8.0, 0.006},
976 { 9.0, 0.0032},
977 { 10.0, 0.0025} };
978
979 sBarkasCorr = new G4PhysicsFreeVector(47, false);
980 for(i=0; i<47; ++i) { sBarkasCorr->PutValues(i, fTable[i][0], fTable[i][1]); }
981
982 static const G4double SK[20] = {1.9477, 1.9232, 1.8996, 1.8550, 1.8137,
983 1.7754, 1.7396, 1.7223, 1.7063, 1.6752,
984 1.6461, 1.6189, 1.5933, 1.5811, 1.5693,
985 1.5467, 1.5254, 1.5053, 1.4863, 1.4772};
986 static const G4double TK[20] = {2.5222, 2.5125, 2.5026, 2.4821, 2.4608,
987 2.4388, 2.4163, 2.4044, 2.3933, 2.3701,
988 2.3466, 2.3229, 2.2992, 2.2872, 2.2753,
989 2.2515, 2.2277, 2.2040, 2.1804, 2.1686};
990
991 static const G4double SL[26] = {15.3343, 13.9389, 12.7909, 11.8343, 11.0283,
992 10.3424, 10.0371, 9.7537, 9.2443, 8.8005,
993 8.4114, 8.0683, 7.9117, 7.7641, 7.4931,
994 7.2506, 7.0327, 6.8362, 6.7452, 6.6584,
995 6.4969, 6.3498, 6.2154, 6.0923, 6.0345, 5.9792};
996 static const G4double TL[26] = {35.0669, 33.4344, 32.0073, 30.7466, 29.6226,
997 28.6128, 28.1449, 27.6991, 26.8674, 26.1061,
998 25.4058, 24.7587, 24.4531, 24.1583, 23.5992,
999 23.0771, 22.5880, 22.1285, 21.9090, 21.6958,
1000 21.2872, 20.9006, 20.5341, 20.1859, 20.0183, 19.8546};
1001
1002 static const G4double bk1[29][11] = {
1003 {0.005, 1.34782E-8, 1.46132E-8, 1.72179E-8, 2.03521E-8, 2.41370E-8, 2.87247E-8, 3.13778E-8, 3.43072E-8, 4.11274E-8, 4.94946E-8},
1004 {0.007, 6.87555E-8, 7.44373E-8, 8.74397E-8, 1.03022E-7, 1.21760E-7, 1.44370E-7, 1.57398E-7, 1.71747E-7, 2.05023E-7, 2.45620E-7},
1005 {0.01, 3.78413E-7, 4.08831E-7, 4.78154E-7, 5.60760E-7, 6.59478E-7, 7.77847E-7, 8.45709E-7, 9.20187E-7, 1.09192E-6, 1.29981E-6},
1006 {0.015, 2.53200E-6, 2.72664E-6, 3.16738E-6, 3.68791E-6, 4.30423E-6, 5.03578E-6, 5.45200E-6, 5.90633E-6, 6.94501E-6, 8.18757E-6},
1007 {0.02, 9.43891E-6, 1.01339E-5, 1.16984E-5, 1.35313E-5, 1.56829E-5, 1.82138E-5, 1.96439E-5, 2.11973E-5, 2.47216E-5, 2.88935E-5},
1008 {0.03, 5.67088E-5, 6.05576E-5, 6.91311E-5, 7.90324E-5, 9.04832E-5, 1.03744E-4, 1.11147E-4, 1.19122E-4, 1.36980E-4, 1.57744E-4},
1009 {0.04, 1.91576E-4, 2.03626E-4, 2.30230E-4, 2.60584E-4, 2.95248E-4, 3.34870E-4, 3.56771E-4, 3.80200E-4, 4.32104E-4, 4.91584E-4},
1010 {0.05, 4.74226E-4, 5.02006E-4, 5.62872E-4, 6.31597E-4, 7.09244E-4, 7.97023E-4, 8.45134E-4, 8.96410E-4, 1.00867E-3, 1.13590E-3},
1011 {0.06, 9.67285E-4, 1.02029E-3, 1.13566E-3, 1.26476E-3, 1.46928E-3, 1.57113E-3, 1.65921E-3, 1.75244E-3, 1.95562E-3, 2.18336E-3},
1012 {0.08, 2.81537E-3, 2.95200E-3, 3.24599E-3, 3.57027E-3, 3.92793E-3, 4.32246E-3, 4.53473E-3, 4.75768E-3, 5.23785E-3, 5.76765E-3},
1013 {0.1, 6.14216E-3, 6.40864E-3, 6.97750E-3, 7.59781E-3, 8.27424E-3, 9.01187E-3, 9.40534E-3, 9.81623E-3, 1.06934E-2, 1.16498E-2},
1014 {0.15, 2.23096E-2, 2.30790E-2, 2.46978E-2, 2.64300E-2, 2.82832E-2, 3.02661E-2, 3.13090E-2, 3.23878E-2, 3.46580E-2, 3.70873E-2},
1015 {0.2, 5.04022E-2, 5.18492E-2, 5.48682E-2, 5.80617E-2, 6.14403E-2, 6.50152E-2, 6.68798E-2, 6.87981E-2, 7.28020E-2, 7.70407E-2},
1016 {0.3, 1.38018E-1, 1.41026E-1, 1.47244E-1, 1.53743E-1, 1.60536E-1, 1.67641E-1, 1.71315E-1, 1.75074E-1, 1.82852E-1, 1.90997E-1},
1017 {0.4, 2.56001E-1, 2.60576E-1, 2.69992E-1, 2.79776E-1, 2.89946E-1, 3.00525E-1, 3.05974E-1, 3.11533E-1, 3.22994E-1, 3.34935E-1},
1018 {0.5, 3.92181E-1, 3.98213E-1, 4.10597E-1, 4.23427E-1, 4.36726E-1, 4.50519E-1, 4.57610E-1, 4.64834E-1, 4.79700E-1, 4.95148E-1},
1019 {0.6, 5.37913E-1, 5.45268E-1, 5.60350E-1, 5.75948E-1, 5.92092E-1, 6.08811E-1, 6.17396E-1, 6.26136E-1, 6.44104E-1, 6.62752E-1},
1020 {0.7, 6.87470E-1, 6.96021E-1, 7.13543E-1, 7.31650E-1, 7.50373E-1, 7.69748E-1, 7.79591E-1, 7.89811E-1, 8.10602E-1, 8.32167E-1},
1021 {0.8, 8.37159E-1, 8.46790E-1, 8.66525E-1, 8.86910E-1, 9.07979E-1, 9.29772E-1, 9.40953E-1, 9.52331E-1, 9.75701E-1, 9.99934E-1},
1022 {1.0, 1.12850, 1.14002, 1.16362, 1.18799, 1.21317, 1.23921, 1.25257, 1.26616, 1.29408, 1.32303},
1023 {1.2, 1.40232, 1.41545, 1.44232, 1.47007, 1.49875, 1.52842, 1.54364, 1.55913, 1.59095, 1.62396},
1024 {1.4, 1.65584, 1.67034, 1.70004, 1.73072, 1.76244, 1.79526, 1.81210, 1.82925, 1.86448, 1.90104},
1025 {1.5, 1.77496, 1.79009, 1.82107, 1.85307, 1.88617, 1.92042, 1.93800, 1.95590, 1.99269, 2.03087},
1026 {1.7, 1.99863, 2.01490, 2.04822, 2.08265, 2.11827, 2.15555, 2.17409, 2.19337, 2.23302, 2.27419},
1027 {2.0, 2.29325, 2.31100, 2.34738, 2.38501, 2.42395, 2.46429, 2.48401, 2.50612, 2.54955, 2.59468},
1028 {3.0, 3.08887, 3.11036, 3.15445, 3.20013, 3.24748, 3.29664, 3.32192, 3.34770, 3.40081, 3.45611},
1029 {5.0, 4.07599, 4.10219, 4.15606, 4.21199, 4.27010, 4.33056, 4.36172, 4.39353, 4.45918, 4.52772},
1030 {7.0, 4.69647, 4.72577, 4.78608, 4.84877, 4.91402, 4.98200, 5.01707, 5.05290, 5.12695, 5.20436},
1031 {10.0, 5.32590, 5.35848, 5.42560, 5.49547, 5.56830, 5.64429, 5.68353, 5.72366, 5.80666, 5.89359}
1032 };
1033
1034 static const G4double bk2[29][11] = {
1035 {0.005, 5.98040E-8, 7.25636E-8, 8.00602E-8, 8.84294E-8, 1.08253E-7, 1.33148E-7, 1.64573E-7, 2.04459E-7, 2.28346E-7, 2.55370E-7},
1036 {0.007, 2.95345E-7, 3.56497E-7, 3.92247E-7, 4.32017E-7, 5.25688E-7, 6.42391E-7, 7.88464E-7, 9.72171E-7, 1.08140E-6, 1.20435E-6},
1037 {0.01, 1.55232E-6, 1.86011E-6, 2.03881E-6, 2.23662E-6, 2.69889E-6, 3.26860E-6, 3.26860E-6, 4.84882E-6, 5.36428E-6, 5.94048E-6},
1038 {0.015, 9.67802E-6, 1.14707E-5, 1.25008E-5, 1.36329E-5, 1.62480E-5, 1.94200E-5, 2.32783E-5, 2.79850E-5, 3.07181E-5, 3.37432E-5},
1039 {0.02, 3.38425E-5, 3.97259E-5, 4.30763E-5, 4.67351E-5, 5.51033E-5, 6.51154E-5, 7.71154E-5, 9.15431E-5, 9.98212E-5, 1.08909E-4},
1040 {0.03, 1.81920E-4, 2.10106E-4, 2.25918E-4, 2.43007E-4, 2.81460E-4, 3.26458E-4, 3.79175E-4, 4.41006E-4, 4.75845E-4, 5.13606E-4},
1041 {0.04, 5.59802E-4, 6.38103E-4, 6.81511E-4, 7.28042E-4, 8.31425E-4, 9.50341E-4, 1.08721E-3, 1.24485E-3, 1.33245E-3, 1.42650E-3},
1042 {0.05, 1.28002E-3, 1.44336E-3, 1.53305E-3, 1.62855E-3, 1.83861E-3, 2.07396E-3, 2.34750E-3, 2.65469E-3, 2.82358E-3, 3.00358E-3},
1043 {0.06, 2.42872E-3, 2.72510E-3, 2.88111E-3, 3.04636E-3, 3.40681E-3, 3.81132E-3, 4.26536E-3, 4.77507E-3, 5.05294E-3, 5.34739E-3},
1044 {0.08, 6.35222E-3, 6.99730E-3, 7.34446E-3, 7.70916E-3, 8.49478E-3, 9.36187E-3, 1.03189E-2, 1.13754E-2, 1.19441E-2, 1.25417E-2},
1045 {0.1, 1.26929E-2, 1.38803E-2, 1.44371E-2, 1.50707E-2, 1.64235E-2, 1.78989E-2, 1.95083E-2, 2.12640E-2, 2.22009E-2, 2.31795E-2},
1046 {0.15, 3.96872E-2, 4.24699E-2, 4.39340E-2, 4.54488E-2, 4.86383E-2, 5.20542E-2, 5.57135E-2, 5.96350E-2, 6.17003E-2, 6.38389E-2},
1047 {0.2, 8.15290E-2, 8.62830E-2, 8.87650E-2, 9.13200E-2, 9.66589E-2, 1.02320E-1, 1.08326E-1, 1.14701E-1, 1.18035E-1, 1.21472E-1},
1048 {0.3, 1.99528E-1, 2.08471E-1, 2.13103E-1, 2.17848E-1, 2.27689E-1, 2.38022E-1, 2.48882E-1, 2.60304E-1, 2.66239E-1, 2.72329E-1},
1049 {0.4, 3.47383E-1, 3.60369E-1, 3.67073E-1, 3.73925E-1, 3.88089E-1, 4.02900E-1, 4.18404E-1, 4.34647E-1, 4.43063E-1, 4.51685E-1},
1050 {0.5, 5.11214E-1, 5.27935E-1, 5.36554E-1, 5.45354E-1, 5.63515E-1, 5.82470E-1, 6.02273E-1, 6.22986E-1, 6.33705E-1, 6.44677E-1},
1051 {0.6, 6.82122E-1, 7.02260E-1, 7.12631E-1, 7.23214E-1, 7.45041E-1, 7.67800E-1, 7.91559E-1, 8.16391E-1, 8.29235E-1, 8.42380E-1},
1052 {0.7, 8.54544E-1, 8.77814E-1, 8.89791E-1, 9.02008E-1, 9.27198E-1, 9.53454E-1, 9.80856E-1, 1.00949, 1.02430, 1.03945},
1053 {0.8, 1.02508, 1.05121, 1.06466, 1.07838, 1.10667, 1.13615, 1.16692, 1.19907, 1.21570, 1.23272},
1054 {1.0, 1.35307, 1.38429, 1.40036, 1.41676, 1.45057, 1.48582, 1.52263, 1.56111, 1.58102, 1.60140},
1055 {1.2, 1.65823, 1.69385, 1.71220, 1.73092, 1.76954, 1.80983, 1.85192, 1.89596, 1.91876, 1.94211},
1056 {1.4, 1.93902, 1.97852, 1.99887, 2.01964, 2.06251, 2.10727, 2.15406, 2.20304, 2.22842, 2.25442},
1057 {1.5, 2.07055, 2.11182, 2.13309, 2.15480, 2.19963, 2.24644, 2.29539, 2.34666, 2.37323, 2.40045},
1058 {1.7, 2.31700, 2.36154, 2.38451, 2.40798, 2.45641, 2.50703, 2.56000, 2.61552, 2.64430, 2.67381},
1059 {2.0, 2.64162, 2.69053, 2.71576, 2.74154, 2.79481, 2.85052, 2.90887, 2.97009, 3.00185, 3.03442},
1060 {3.0, 3.51376, 3.57394, 3.60503, 3.63684, 3.70268, 3.77170, 3.84418, 3.92040, 3.96003, 4.00073},
1061 {5.0, 4.59935, 4.67433, 4.71316, 4.75293, 4.83543, 4.92217, 5.01353, 5.10992, 5.16014, 5.21181},
1062 {7.0, 5.28542, 5.37040, 5.41445, 5.45962, 5.55344, 5.65226, 5.79496, 5.90517, 5.96269, 6.02191},
1063 {10.0, 5.98474, 6.08046, 6.13015, 6.18112, 6.28715, 6.39903, 6.51728, 6.64249, 6.70792, 6.77535}
1064 };
1065
1066 static const G4double bls1[28][10] = {
1067 {0.005, 2.4111E-4, 2.5612E-4, 2.7202E-4, 3.0658E-4, 3.4511E-4, 3.8795E-4, 4.3542E-4, 4.6100E-4, 4.8786E-4},
1068 {0.007, 6.3947E-4, 6.7058E-4, 7.0295E-4, 7.7167E-4, 8.4592E-4, 9.2605E-4, 1.0125E-3, 1.0583E-3, 1.1058E-3},
1069 {0.01, 1.5469E-3, 1.6036E-3, 1.6622E-3, 1.7856E-3, 1.9181E-3, 2.1615E-3, 2.3178E-3, 2.4019E-3, 2.4904E-3},
1070 {0.015, 3.7221E-3, 3.8404E-3, 3.9650E-3, 4.2354E-3, 4.5396E-3, 4.8853E-3, 5.2820E-3, 5.5031E-3, 5.7414E-3},
1071 {0.02, 6.9449E-3, 7.1910E-3, 7.4535E-3, 8.0336E-3, 8.6984E-3, 9.4638E-3, 1.0348E-2, 1.0841E-2, 1.1372E-2},
1072 {0.03, 1.7411E-2, 1.8180E-2, 1.8997E-2, 2.0784E-2, 2.2797E-2, 2.5066E-2, 2.7622E-2, 2.9020E-2, 3.0503E-2},
1073 {0.04, 3.8474E-2, 4.0056E-2, 4.1718E-2, 4.5300E-2, 4.9254E-2, 5.3619E-2, 5.8436E-2, 6.1028E-2, 6.3752E-2},
1074 {0.05, 6.7371E-2, 6.9928E-2, 7.2596E-2, 7.8282E-2, 8.4470E-2, 9.1206E-2, 9.8538E-2, 1.0244E-1, 1.0652E-1},
1075 {0.06, 1.0418E-1, 1.0778E-1, 1.1152E-1, 1.1943E-1, 1.2796E-1, 1.3715E-1, 1.4706E-1, 1.5231E-1, 1.5776E-1},
1076 {0.08, 1.9647E-1, 2.0217E-1, 2.0805E-1, 2.2038E-1, 2.3351E-1, 2.4751E-1, 2.6244E-1, 2.7027E-1, 2.7837E-1},
1077 {0.1, 3.0594E-1, 3.1361E-1, 3.2150E-1, 3.3796E-1, 3.5537E-1, 3.7381E-1, 3.9336E-1, 4.0357E-1, 4.1410E-1},
1078 {0.15, 6.1411E-1, 6.2597E-1, 6.3811E-1, 6.6330E-1, 6.8974E-1, 7.1753E-1, 7.4678E-1, 7.6199E-1, 7.7761E-1},
1079 {0.2, 9.3100E-1, 9.5614E-1, 9.7162E-1, 1.0037, 1.0372, 1.0723, 1.1092, 1.1284, 1.1480},
1080 {0.3, 1.5172, 1.5372, 1.5576, 1.5998, 1.6438, 1.6899, 1.7382, 1.7632, 1.7889},
1081 {0.4, 2.0173, 2.0408, 2.0647, 2.1142, 2.1659, 2.2199, 2.2765, 2.3059, 2.3360},
1082 {0.5, 2.3932, 2.4193, 2.4460, 2.5011, 2.5587, 2.6190, 2.6821, 2.7148, 2.7484},
1083 {0.6, 2.7091, 2.7374, 2.7663, 2.8260, 2.8884, 2.9538, 3.0222, 3.0577, 3.0941},
1084 {0.7, 2.9742, 3.0044, 3.0352, 3.0988, 3.1652, 3.2349, 3.3079, 3.3457, 3.3845},
1085 {0.8, 3.2222, 3.2539, 3.2863, 3.3532, 3.4232, 3.4965, 3.5734, 3.6133, 3.6542},
1086 {1.0, 3.6690, 3.7033, 3.7384, 3.8108, 3.8866, 3.9661, 4.0495, 4.0928, 4.1371},
1087 {1.2, 3.9819, 4.0183, 4.0555, 4.1324, 4.2130, 4.2974, 4.3861, 4.4321, 4.4794},
1088 {1.4, 4.2745, 4.3127, 4.3517, 4.4324, 4.5170, 4.6056, 4.6988, 4.7471, 4.7968},
1089 {1.5, 4.4047, 4.4436, 4.4834, 4.5658, 4.6521, 4.7426, 4.8378, 4.8872, 4.9379},
1090 {1.7, 4.6383, 4.6787, 4.7200, 4.8054, 4.8949, 4.9888, 5.0876, 5.1388, 5.1915},
1091 {2.0, 4.9369, 4.9791, 5.0223, 5.1116, 5.2053, 5.3036, 5.4070, 5.4607, 5.5159},
1092 {3.0, 5.6514, 5.6981, 5.7459, 5.8450, 5.9489, 6.0581, 6.1730, 6.2328, 6.2943},
1093 {5.0, 6.4665, 6.5189, 6.5724, 6.6835, 6.8003, 6.9231, 7.0525, 7.1199, 7.1892},
1094 {7.0, 6.8634, 6.9194, 6.9767, 7.0957, 7.2208, 7.3526, 7.4915, 7.5639, 7.6384}
1095 };
1096
1097 static const G4double bls2[28][10] = {
1098 {0.005, 5.4561E-4, 6.0905E-4, 6.7863E-4, 7.5494E-4, 7.9587E-4, 8.3883E-4, 9.3160E-4, 1.0352E-3, 1.1529E-3},
1099 {0.007, 1.2068E-3, 1.3170E-3, 1.4377E-3, 1.5719E-3, 1.6451E-3, 1.7231E-3, 1.8969E-3, 2.1009E-3, 2.3459E-3},
1100 {0.01, 2.6832E-3, 2.9017E-3, 3.1534E-3, 3.4479E-3, 3.6149E-3, 3.7976E-3, 4.2187E-3, 4.7320E-3, 5.3636E-3},
1101 {0.015, 6.2775E-3, 6.9077E-3, 7.6525E-3, 8.5370E-2, 9.0407E-3, 9.5910E-3, 1.0850E-2, 1.2358E-2, 1.4165E-2},
1102 {0.02, 1.2561E-2, 1.3943E-2, 1.5553E-2, 1.7327E-2, 1.8478E-2, 1.9612E-2, 2.2160E-2, 2.5130E-2, 2.8594E-2},
1103 {0.03, 3.3750E-2, 3.7470E-2, 4.1528E-2, 4.6170E-2, 4.8708E-2, 5.1401E-2, 5.7297E-2, 6.3943E-2, 7.1441E-2},
1104 {0.04, 6.9619E-2, 7.6098E-2, 8.3249E-2, 9.1150E-2, 9.5406E-2, 9.9881E-2, 1.0954E-1, 1.2023E-1, 1.3208E-1},
1105 {0.05, 1.1522E-1, 1.2470E-1, 1.3504E-1, 1.4632E-1, 1.5234E-1, 1.5864E-1, 1.7211E-1, 1.8686E-1, 2.0304E-1},
1106 {0.06, 1.6931E-1, 1.8179E-1, 1.9530E-1, 2.0991E-1, 2.1767E-1, 2.2576E-1, 2.4295E-1, 2.6165E-1, 2.8201E-1},
1107 {0.08, 2.9540E-1, 3.1361E-1, 3.3312E-1, 3.5403E-1, 3.6506E-1, 3.7650E-1, 4.0067E-1, 4.2673E-1, 4.5488E-1},
1108 {0.1, 4.3613E-1, 4.5956E-1, 4.852E-1, 5.1115E-1, 5.2514E-1, 5.3961E-1, 5.7008E-1, 6.0277E-1, 6.3793E-1},
1109 {0.15, 8.1014E-1, 8.4453E-1, 8.8093E-1, 9.1954E-1, 9.3973E-1, 9.6056E-1, 1.0043, 1.0509, 1.1008},
1110 {0.2, 1.1888, 1.2319, 1.2774, 1.3255, 1.3506, 1.3765, 1.4308, 1.4886, 1.5504},
1111 {0.3, 1.8422, 1.8983, 1.9575, 2.0201, 2.0528, 2.0864, 2.1569, 2.2319, 2.3120},
1112 {0.4, 2.3984, 2.4642, 2.5336, 2.6070, 2.6452, 2.6847, 2.7674, 2.8554, 2.9494},
1113 {0.5, 2.8181, 2.8915, 2.9690, 3.0509, 3.0937, 3.1378, 3.2301, 3.3285, 3.4337},
1114 {0.6, 3.1698, 3.2494, 3.3336, 3.4226, 3.4691, 3.5171, 3.6175, 3.7246, 3.8391},
1115 {0.7, 3.4652, 3.5502, 3.6400, 3.7351, 3.7848, 3.8360, 3.9433, 4.0578, 4.1804},
1116 {0.8, 3.7392, 3.8289, 3.9236, 4.0239, 4.0764, 4.1304, 4.2438, 4.3648, 4.4944},
1117 {1.0, 4.2295, 4.3269, 4.4299, 4.5391, 4.5962, 4.6551, 4.7786, 4.9106, 5.0520},
1118 {1.2, 4.5777, 4.6814, 4.7912, 4.9076, 4.9685, 5.0314, 5.1633, 5.3043, 5.4555},
1119 {1.4, 4.9001, 5.0092, 5.1247, 5.2473, 5.3114, 5.3776, 5.5166, 5.6653, 5.8249},
1120 {1.5, 5.0434, 5.1550, 5.2731, 5.3984, 5.4640, 5.5317, 5.6739, 5.8260, 5.9893},
1121 {1.7, 5.3011, 5.4170, 5.5398, 5.6701, 5.7373, 5.8088, 5.9568, 6.1152, 6.2853},
1122 {2.0, 5.6308, 5.7523, 5.8811, 6.0174, 6.0896, 6.1636, 6.3192, 6.4857, 6.6647},
1123 {3.0, 6.4224, 6.5580, 6.7019, 6.8549, 6.9351, 7.0180, 7.1925, 7.3795, 7.5808},
1124 {5.0, 7.3338, 7.4872, 7.6500, 7.8235, 7.9146, 8.0087, 8.2071, 8.4200, 8.6496},
1125 {7.0, 7.7938, 7.9588, 8.1342, 8.3211, 8.4193, 8.5209, 8.7350, 8.9651, 9.2133}
1126 };
1127
1128 static const G4double bls3[28][9] = {
1129 {0.005, 1.2895E-3, 1.3670E-3, 1.4524E-3, 1.6524E-3, 1.9078E-3, 2.2414E-3, 2.6889E-3, 3.3006E-3},
1130 {0.007, 2.6467E-3, 2.8242E-3, 3.0238E-3, 3.5045E-3, 4.1260E-3, 4.9376E-3, 6.0050E-3, 7.4152E-3},
1131 {0.01, 6.1472E-3, 6.6086E-3, 7.1246E-3, 8.3491E-3, 9.8871E-3, 1.1822E-2, 1.4261E-2, 1.7335E-2},
1132 {0.015, 1.63316E-2, 1.7572E-2, 1.8932E-2, 2.2053E-2, 2.5803E-2, 3.0308E-2, 3.5728E-2, 4.2258E-2},
1133 {0.02, 3.2634E-2, 3.4900E-2, 3.7348E-2, 4.2850E-2, 4.9278E-2, 5.6798E-2, 6.5610E-2, 7.5963E-2},
1134 {0.03, 7.9907E-2, 8.4544E-2, 8.9486E-2, 1.0032E-1, 1.1260E-1, 1.2656E-1, 1.4248E-1, 1.6071E-1},
1135 {0.04, 1.4523E-1, 1.5237E-1, 1.5985E-1, 1.7614E-1, 1.9434E-1, 2.1473E-1, 2.3766E-1, 2.6357E-1},
1136 {0.05, 2.2082E-1, 2.3036E-1, 2.4038E-1, 2.6199E-1, 2.8590E-1, 3.1248E-1, 3.4212E-1, 3.7536E-1},
1137 {0.06, 3.0423E-1, 3.1611E-1, 3.2854E-1, 3.5522E-1, 3.8459E-1, 4.1704E-1, 4.5306E-1, 4.9326E-1},
1138 {0.08, 4.8536E-1, 5.0156E-1, 5.1846E-1, 5.5453E-1, 5.9397E-1, 6.3728E-1, 6.8507E-1, 7.3810E-1},
1139 {0.1, 6.7586E-1, 6.9596E-1, 7.1688E-1, 7.6141E-1, 8.0992E-1, 8.6301E-1, 9.2142E-1, 9.8604E-1},
1140 {0.15, 1.1544, 1.1828, 1.2122, 1.2746, 1.3424, 1.4163, 1.4974, 1.5868},
1141 {0.2, 1.6167, 1.6517, 1.6880, 1.7650, 1.8484, 1.9394, 2.0390, 2.1489},
1142 {0.3, 2.3979, 2.4432, 2.4902, 2.5899, 2.6980, 2.8159, 2.9451, 3.0876},
1143 {0.4, 3.0502, 3.1034, 3.1586, 3.2758, 3.4030, 3.5416, 3.6938, 3.8620},
1144 {0.5, 3.5466, 3.6062, 3.6681, 3.7994, 3.9421, 4.0978, 4.2688, 4.4580},
1145 {0.6, 3.9620, 4.0270, 4.0945, 4.2378, 4.3935, 4.5636, 4.7506, 4.9576},
1146 {0.7, 4.3020, 4.3715, 4.4438, 4.5974, 4.7644, 4.9470, 5.1478, 5.3703},
1147 {0.8, 4.6336, 4.7072, 4.7837, 4.9463, 5.1233, 5.3169, 5.5300, 5.7661},
1148 {1.0, 5.2041, 5.2845, 5.3682, 5.5462, 5.7400, 5.9523, 6.1863, 6.4458},
1149 {1.2, 5.6182, 5.7044, 5.7940, 5.9848, 6.1927, 6.4206, 6.6719, 6.9510},
1150 {1.4, 5.9967, 6.0876, 6.1823, 6.3839, 6.6038, 6.8451, 7.1113, 7.4071},
1151 {1.5, 6.1652, 6.2583, 6.3553, 6.5618, 6.7871, 7.0344, 7.3073, 7.6107},
1152 {1.7, 6.4686, 6.5657, 6.6668, 6.8823, 7.1175, 7.3757, 7.6609, 7.9782},
1153 {2.0, 6.8577, 6.9600, 7.0666, 7.2937, 7.5418, 7.8143, 8.1156, 8.4510},
1154 {3.0, 7.7981, 7.9134, 8.0336, 8.2901, 8.5708, 8.8796, 9.2214, 9.6027},
1155 {5.0, 8.8978, 9.0297, 9.1673, 9.4612, 9.7834, 10.1384, 10.5323, 10.9722},
1156 {7.0, 9.4819, 9.6248, 9.7739, 10.0926, 10.4423, 10.8282, 11.2565, 11.7356}
1157 };
1158
1159 static const G4double bll1[28][10] = {
1160 {0.005, 3.6324E-5, 4.0609E-5, 4.5430E-5, 5.6969E-5, 7.1625E-5, 9.0279E-5, 1.1407E-4, 1.2834E-4, 1.4447E-4},
1161 {0.007, 1.8110E-4, 2.0001E-4, 2.2099E-4, 2.7006E-4, 3.3049E-4, 4.0498E-4, 4.9688E-4, 5.5061E-4, 6.1032E-4},
1162 {0.01, 8.6524E-4, 9.4223E-4, 1.0262E-3, 1.2178E-3, 1.4459E-3, 1.7174E-3, 2.0405E-3, 2.2245E-3, 2.4252E-3},
1163 {0.015, 4.2293E-3, 4.5314E-3, 4.8551E-3, 5.5731E-3, 6.3968E-3, 7.3414E-3, 8.4242E-3, 9.0236E-3, 9.6652E-3},
1164 {0.02, 1.1485E-2, 1.2172E-2, 1.2900E-2, 1.4486E-2, 1.6264E-2, 1.8256E-2, 2.0487E-2, 2.1702E-2, 2.2989E-2},
1165 {0.03, 3.9471E-2, 4.1270E-2, 4.3149E-2, 4.7163E-2, 5.1543E-2, 5.6423E-2, 6.1540E-2, 6.4326E-2, 6.7237E-2},
1166 {0.04, 8.4454E-2, 8.7599E-2, 9.0860E-2, 9.7747E-2, 1.0516E-1, 1.1313E-1, 1.2171E-1, 1.2625E-1, 1.3096E-1},
1167 {0.05, 1.4339E-1, 1.4795E-1, 1.5266E-1, 1.6253E-1, 1.7306E-1, 1.8430E-1, 1.9630E-1, 2.0261E-1, 2.0924E-1},
1168 {0.06, 2.1304E-1, 2.1899E-1, 2.2512E-1, 2.3794E-1, 2.5153E-1, 2.6596E-1, 2.8130E-1, 2.8934E-1, 2.9763E-1},
1169 {0.08, 3.7382E-1, 3.8241E-1, 3.9122E-1, 4.0955E-1, 4.2888E-1, 4.4928E-1, 4.7086E-1, 4.8212E-1, 4.9371E-1},
1170 {0.1, 5.5056E-1, 5.6151E-1, 5.7273E-1, 5.9601E-1, 6.2049E-1, 6.4627E-1, 6.7346E-1, 6.8762E-1, 7.0218E-1},
1171 {0.15, 1.0066, 1.0224, 1.0386, 1.0721, 1.1073, 1.1443, 1.1832, 1.2035, 1.2243},
1172 {0.2, 1.4376, 1.4572, 1.4773, 1.5188, 1.5624, 1.6083, 1.6566, 1.6817, 1.7076},
1173 {0.3, 2.1712, 2.1964, 2.2223, 2.2758, 2.3322, 2.3915, 2.4542, 2.4869, 2.5205},
1174 {0.4, 2.7500, 2.7793, 2.8094, 2.8719, 2.9377, 3.0072, 3.0807, 3.1192, 3.1587},
1175 {0.5, 3.2033, 3.2359, 3.2693, 3.3389, 3.4122, 3.4898, 3.5721, 3.6151, 3.6595},
1176 {0.6, 3.6038, 3.6391, 3.6753, 3.7506, 3.8303, 3.9146, 4.0042, 4.0511, 4.0995},
1177 {0.7, 3.9106, 3.9482, 3.9867, 4.0670, 4.1520, 4.2421, 4.3380, 4.3882, 4.4401},
1178 {0.8, 4.1790, 4.2185, 4.2590, 4.3437, 4.4333, 4.5285, 4.6298, 4.6830, 4.7380},
1179 {1.0, 4.6344, 4.6772, 4.7212, 4.8131, 4.9106, 5.0144, 5.1250, 5.1831, 5.2432},
1180 {1.2, 4.9787, 5.0242, 5.0711, 5.1689, 5.2729, 5.3837, 5.5050, 5.5642, 5.6287},
1181 {1.4, 5.2688, 5.3166, 5.3658, 5.4687, 5.5782, 5.6950, 5.8198, 5.8855, 5.9554},
1182 {1.5, 5.3966, 5.4454, 5.4957, 5.6009, 5.7128, 5.8323, 5.9601, 6.0274, 6.0972},
1183 {1.7, 5.6255, 5.6762, 5.7284, 5.8377, 5.9541, 6.0785, 6.2116, 6.2818, 6.3546},
1184 {2.0, 5.9170, 5.9701, 6.0248, 6.1395, 6.2618, 6.3925, 6.5327, 6.6066, 6.6833},
1185 {3.0, 6.6210, 6.6801, 6.7411, 6.8692, 7.0062, 7.1529, 7.3107, 7.3941, 7.4807},
1186 {5.0, 7.4620, 7.5288, 7.5977, 7.7428, 7.8982, 8.0653, 8.2454, 8.3409, 8.4402},
1187 {7.0, 7.7362, 7.8079, 7.8821, 8.0383, 8.2061, 8.3866, 8.5816, 8.6850, 8.7927}
1188 };
1189
1190 static const G4double bll2[28][10] = {
1191 {0.005, 1.8339E-4, 2.3330E-4, 2.9738E-4, 3.7977E-4, 4.2945E-4, 4.8582E-4, 6.2244E-4, 7.9858E-4, 1.0258E-3},
1192 {0.007, 7.5042E-4, 9.2355E-4, 1.1375E-3, 1.4021E-3, 1.5570E-3, 1.7292E-3, 2.1335E-3, 2.6335E-3, 3.2515E-3},
1193 {0.01, 2.8829E-3, 3.4275E-3, 4.0758E-3, 4.8457E-3, 5.2839E-3, 5.7617E-3, 6.8504E-3, 8.1442E-3, 9.6816E-3},
1194 {0.015, 1.1087E-2, 1.2716E-2, 1.4581E-2, 1.6717E-2, 1.7898E-2, 1.9163E-2, 2.1964E-2, 2.5173E-2, 2.8851E-2},
1195 {0.02, 2.5786E-2, 2.8922E-2, 3.2435E-2, 3.6371E-2, 3.8514E-2, 4.0784E-2, 4.5733E-2, 5.1288E-2, 5.7531E-2},
1196 {0.03, 7.3461E-2, 8.0264E-2, 8.7705E-2, 9.5852E-2, 1.0021E-1, 1.0478E-1, 1.1458E-1, 1.2535E-1, 1.3721E-1},
1197 {0.04, 1.4094E-1, 1.5172E-1, 1.6336E-1, 1.7596E-1, 1.8265E-1, 1.8962E-1, 2.0445E-1, 2.2058E-1, 2.3818E-1},
1198 {0.05, 2.2289E-1, 2.3762E-1, 2.5344E-1, 2.7045E-1, 2.7944E-1, 2.8877E-1, 3.0855E-1, 3.2995E-1, 3.5318E-1},
1199 {0.06, 3.1503E-1, 3.3361E-1, 3.5346E-1, 3.7473E-1, 3.8594E-1, 3.9756E-1, 4.2212E-1, 4.4861E-1, 4.7727E-1},
1200 {0.08, 5.1793E-1, 5.4368E-1, 5.7109E-1, 6.0032E-1, 6.1569E-1, 6.3159E-1, 6.6512E-1, 7.0119E-1, 7.4012E-1},
1201 {0.1, 7.3258E-1, 7.6481E-1, 7.9907E-1, 8.3556E-1, 8.5472E-1, 8.7454E-1, 9.1630E-1, 9.6119E-1, 1.0096},
1202 {0.15, 1.2677, 1.3137, 1.3626, 1.4147, 1.4421, 1.4703, 1.5300, 1.5942, 1.6636},
1203 {0.2, 1.7615, 1.8188, 1.8797, 1.9446, 1.9788, 2.0142, 2.0889, 2.1694, 2.2567},
1204 {0.3, 2.5909, 2.6658, 2.7457, 2.8312, 2.8762, 2.9231, 3.0222, 3.1295, 3.2463},
1205 {0.4, 3.2417, 3.3302, 3.4249, 3.5266, 3.5803, 3.6361, 3.7546, 3.8835, 4.0242},
1206 {0.5, 3.7527, 3.8523, 3.9591, 4.0741, 4.1350, 4.1983, 4.3330, 4.4799, 4.6408},
1207 {0.6, 4.2013, 4.3103, 4.4274, 4.5537, 4.6206, 4.6904, 4.8390, 5.0013, 5.1796},
1208 {0.7, 4.5493, 4.6664, 4.7925, 4.9286, 5.0009, 5.0762, 5.2370, 5.4129, 5.6066},
1209 {0.8, 4.8537, 4.9780, 5.1119, 5.2568, 5.3338, 5.4141, 5.5857, 5.7738, 5.9811},
1210 {1.0, 5.3701, 5.5066, 5.6540, 5.8138, 5.8989, 5.9878, 6.1780, 6.3870, 6.6179},
1211 {1.2, 5.7648, 5.9114, 6.0701, 6.2424, 6.3343, 6.4303, 6.6361, 6.8626, 7.1137},
1212 {1.4, 6.0976, 6.2530, 6.4214, 6.6044, 6.7021, 6.8043, 7.0237, 7.2655, 7.5338},
1213 {1.5, 6.2447, 6.4041, 6.5768, 6.7647, 6.8650, 6.9700, 7.1954, 7.4442, 7.7203},
1214 {1.7, 6.5087, 6.6752, 6.8558, 7.0526, 7.1578, 7.2679, 7.5045, 7.7660, 8.0565},
1215 {2.0, 6.8458, 7.0218, 7.2129, 7.4213, 7.5328, 7.6496, 7.9010, 8.1791, 8.4886},
1216 {3.0, 7.6647, 7.8644, 8.0819, 8.3189, 8.4475, 8.5814, 8.8702, 9.1908, 9.5488},
1217 {5.0, 8.6515, 8.8816, 9.1330, 9.4090, 9.5572, 9.7132, 10.0504, 10.4259, 10.8466},
1218 {7.0, 9.0221, 9.2724, 9.5464, 9.8477, 10.0099, 10.1805, 10.5499, 10.9622, 11.4250}
1219 };
1220
1221 static const G4double bll3[28][9] = {
1222 {0.005, 1.3190E-3, 1.4961E-3, 1.6974E-3, 2.1858E-3, 2.8163E-3, 3.6302E-3, 4.6814E-3, 6.0395E-3},
1223 {0.007, 4.0158E-3, 4.4623E-3, 4.9592E-3, 6.1257E-3, 7.5675E-3, 9.3502E-3, 1.1556E-2, 1.4290E-2},
1224 {0.01, 1.1509E-2, 1.2548E-2, 1.3681E-2, 1.6263E-2, 1.9336E-2, 2.2999E-2, 2.7370E-2, 3.2603E-2},
1225 {0.015, 3.3070E-2, 3.5408E-2, 3.7914E-2, 4.3483E-2, 4.9898E-2, 5.7304E-2, 6.5884E-2, 7.5861E-2},
1226 {0.02, 6.4555E-2, 6.8394E-2, 7.2472E-2, 8.1413E-2, 9.1539E-2, 1.0304E-1, 1.1617E-1, 1.3121E-1},
1227 {0.03, 1.5030E-1, 1.5101E-1, 1.5844E-1, 1.7451E-1, 1.9244E-1, 2.1244E-1, 2.3496E-1, 2.6044E-1},
1228 {0.04, 2.5743E-1, 2.6774E-1, 2.7855E-1, 3.0180E-1, 3.2751E-1, 3.5608E-1, 3.8803E-1, 4.2401E-1},
1229 {0.05, 3.7846E-1, 3.9195E-1, 4.0607E-1, 4.3635E-1, 4.6973E-1, 5.0672E-1, 5.4798E-1, 5.9436E-1},
1230 {0.06, 5.0839E-1, 5.2497E-1, 5.4230E-1, 5.7943E-1, 6.2028E-1, 6.6549E-1, 7.1589E-1, 7.7252E-1},
1231 {0.08, 7.8230E-1, 8.0474E-1, 8.2818E-1, 8.7836E-1, 9.3355E-1, 9.9462E-1, 1.0627, 1.1394},
1232 {0.1, 1.0621, 1.0900, 1.1192, 1.1816, 1.2503, 1.3265, 1.4116, 1.5076},
1233 {0.15, 1.7389, 1.7790, 1.8210, 1.9112, 2.0108, 2.1217, 2.2462, 2.3876},
1234 {0.2, 2.3516, 2.4024, 2.4556, 2.5701, 2.6971, 2.8391, 2.9994, 3.1822},
1235 {0.3, 3.3741, 3.4427, 3.5148, 3.6706, 3.8445, 4.0404, 4.2631, 4.5193},
1236 {0.4, 4.1788, 4.2620, 4.3496, 4.5398, 4.7530, 4.9944, 5.2703, 5.5895},
1237 {0.5, 4.8180, 4.9137, 5.0146, 5.2341, 5.4811, 5.7618, 6.0840, 6.4583},
1238 {0.6, 5.3765, 5.4830, 5.5954, 5.8406, 6.1173, 6.4326, 6.7958, 7.2191},
1239 {0.7, 5.8208, 5.9369, 6.0596, 6.3275, 6.6306, 6.9769, 7.3767, 7.8440},
1240 {0.8, 6.2109, 6.3355, 6.4674, 6.7558, 7.0827, 7.4570, 7.8900, 8.3972},
1241 {1.0, 6.8747, 7.0142, 7.1621, 7.4861, 7.8546, 8.2778, 8.7690, 9.3464},
1242 {1.2, 7.3933, 7.5454, 7.7068, 8.0612, 8.4652, 8.9302, 9.4713, 10.1090},
1243 {1.4, 7.8331, 7.9967, 8.1694, 8.5502, 8.9851, 9.4866, 10.0713, 10.7619},
1244 {1.5, 8.0286, 8.1967, 8.3753, 8.7681, 9.2181, 9.7352, 10.3399, 11.0546},
1245 {1.7, 8.3813, 8.5585, 8.7469, 9.1618, 9.6367, 10.1856, 10.8270, 11.5863},
1246 {2.0, 8.8352, 9.0245, 9.2260, 9.6701, 10.1793, 10.7688, 11.4590, 12.2775},
1247 {3.0, 9.9511, 10.1714, 10.4062, 10.9254, 11.5229, 12.2172, 13.0332, 14.0048},
1248 {5.0, 11.3211, 11.5818, 11.8601, 12.4771, 13.1898, 14.0213, 15.0024, 16.1752},
1249 {7.0, 11.9480, 12.2357, 12.5432, 13.2260, 14.0164, 14.9404, 16.0330, 17.3420}
1250 };
1251
1252 G4double b, bs;
1253 for(i=0; i<nEtaK; ++i) {
1254
1255 G4double et = Eta[i];
1256 G4double loget = G4Log(et);
1257
1258 for(j=0; j<nK; ++j) {
1259
1260 if(j < 10) { b = bk2[i][10-j]; }
1261 else { b = bk1[i][20-j]; }
1262
1263 CK[j][i] = SK[j]*loget + TK[j] - b;
1264
1265 if(i == nEtaK-1) {
1266 ZK[j] = et*(et*et*CK[j][i] - et*UK[j] - VK[j]);
1267 //G4cout << "i= " << i << " j= " << j
1268 // << " CK[j][i]= " << CK[j][i]
1269 // << " ZK[j]= " << ZK[j] << " b= " << b << G4endl;
1270 }
1271 }
1272 //G4cout << G4endl;
1273 if(i < nEtaL) {
1274 //G4cout << " LShell:" <<G4endl;
1275 for(j=0; j<nL; ++j) {
1276
1277 if(j < 8) {
1278 bs = bls3[i][8-j];
1279 b = bll3[i][8-j];
1280 } else if(j < 17) {
1281 bs = bls2[i][17-j];
1282 b = bll2[i][17-j];
1283 } else {
1284 bs = bls1[i][26-j];
1285 b = bll1[i][26-j];
1286 }
1287 G4double c = SL[j]*loget + TL[j];
1288 CL[j][i] = c - bs - 3.0*b;
1289 if(i == nEtaL-1) {
1290 VL[j] = et*(et*CL[j][i] - UL[j]);
1291 //G4cout << "i= " << i << " j= " << j
1292 // << " CL[j][i]= " << CL[j][i]
1293 // << " VL[j]= " << VL[j] << " b= " << b << " bs= " << bs
1294 // << " et= " << et << G4endl;
1295 //" UL= " << UL[j] << " TL= " << TL[j] << " SL= " << SL[j] <<G4endl;
1296 }
1297 }
1298 //G4cout << G4endl;
1299 }
1300 }
1301
1302 static const G4double xzk[34] = { 11.7711,
1303 13.3669, 15.5762, 17.1715, 18.7667, 20.8523, 23.0606, 24.901, 26.9861, 29.4394, 31.77,
1304 34.3457, 37.4119, 40.3555, 42.3177, 44.7705, 47.2234, 50.78, 53.8458, 56.4214, 58.3834,
1305 60.9586, 63.6567, 66.5998, 68.807, 71.8728, 74.5706, 77.3911, 81.8056, 85.7297, 89.8988,
1306 93.4549, 96.2753, 99.709};
1307 static const G4double yzk[34] = { 0.70663,
1308 0.72033, 0.73651, 0.74647, 0.75518, 0.76388, 0.77258, 0.78129, 0.78625, 0.7937, 0.79991,
1309 0.80611, 0.8123, 0.8185, 0.82097, 0.82467, 0.82838, 0.83457, 0.83702, 0.84198, 0.8432,
1310 0.84565, 0.84936, 0.85181, 0.85303, 0.85548, 0.85794, 0.8604, 0.86283, 0.86527, 0.86646,
1311 0.86891, 0.87011, 0.87381};
1312
1313 static const G4double xzl[36] = { 15.5102,
1314 16.7347, 17.9592, 19.551, 21.0204, 22.6122, 24.9388, 27.3878, 29.5918, 31.3061, 32.898,
1315 34.4898, 36.2041, 38.4082, 40.3674, 42.5714, 44.898, 47.4694, 49.9184, 52.7347, 55.9184,
1316 59.3469, 61.9184, 64.6122, 67.4286, 71.4694, 75.2653, 78.3265, 81.2653, 85.551, 88.7347,
1317 91.551, 94.2449, 96.449, 98.4082, 99.7551};
1318 static const G4double yzl[36] = { 0.29875,
1319 0.31746, 0.33368, 0.35239, 0.36985, 0.38732, 0.41102, 0.43472, 0.45343, 0.4659, 0.47713,
1320 0.4896, 0.50083, 0.51331, 0.52328, 0.53077, 0.54075, 0.54823, 0.55572, 0.56445, 0.57193,
1321 0.58191, 0.5869, 0.59189, 0.60062, 0.60686, 0.61435, 0.61809, 0.62183, 0.62931, 0.6343,
1322 0.6368, 0.64054, 0.64304, 0.64428, 0.64678};
1323
1324 sThetaK = new G4PhysicsFreeVector(34, false);
1325 for(i=0; i<34; ++i) { sThetaK->PutValues(i, xzk[i], yzk[i]); }
1326 sThetaL = new G4PhysicsFreeVector(36, false);
1327 for(i=0; i<36; ++i) { sThetaL->PutValues(i, xzl[i], yzl[i]); }
1328}
1329
1330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1331
static const G4double e1[44]
const G4double alpha2
const G4double inveplus
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static const G4double * SL[nLA]
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
static G4double VL[26]
G4double LShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4VEmModel * ionLEModel
static G4double CK[20][29]
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
std::vector< const G4ParticleDefinition * > ionList
G4VEmModel * ionHEModel
const G4Material * material
static G4double ZK[20]
G4double MottCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ParticleDefinition * curParticle
G4double KShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4Material * curMaterial
G4PhysicsVector * curVector
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
void BuildCorrectionVector()
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double Value(G4double xv, G4double x1, G4double x2, G4double y1, G4double y2) const
std::vector< G4String > materialName
G4double ComputeIonCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4IonTable * ionTable
G4double Value2(G4double xv, G4double yv, G4double x1, G4double x2, G4double y1, G4double y2, G4double z11, G4double z21, G4double z12, G4double z22) const
G4double SpinCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static const G4double ZD[11]
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double KShell(G4double theta, G4double eta)
std::map< G4int, std::vector< G4double > > thcorr
static const G4double Eta[29]
std::vector< G4PhysicsVector * > stopData
static const G4double UL[26]
static const G4double UK[20]
G4int Index(G4double x, const G4double *y, G4int n) const
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ElementVector * theElementVector
G4double Bethe(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ParticleDefinition * particle
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
std::vector< G4int > Zion
G4double ShellCorrectionSTD(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4double CL[26][28]
const G4double * atomDensity
static G4PhysicsFreeVector * sBarkasCorr
void SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double LShell(G4double theta, G4double eta)
static G4PhysicsFreeVector * sThetaK
G4ionEffectiveCharge effCharge
G4double DensityCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
std::vector< const G4Material * > materialList
G4EmCorrections(G4int verb)
static const G4double VK[20]
static G4PhysicsFreeVector * sThetaL
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
std::vector< G4int > Aion
std::vector< const G4Material * > currmat
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4double GetMdensity() const
G4double GetX1density() const
G4double * GetShellCorrectionVector() const
G4double GetTaul() const
G4double GetX0density() const
G4double GetCdensity() const
G4double GetMeanExcitationEnergy() const
G4double GetAdensity() const
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:212
G4double GetElectronDensity() const
Definition: G4Material.hh:213
const G4String & GetName() const
Definition: G4Material.hh:173
G4int GetAtomicNumber() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
static constexpr double eplus
static constexpr double electron_mass_c2
static constexpr double proton_mass_c2
static constexpr double GeV
static constexpr double fine_structure_const
static constexpr double keV
static constexpr double MeV
static constexpr double pi
Definition: SystemOfUnits.h:55
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
float electron_mass_c2
Definition: hepunit.py:273
int twopi_mc2_rcl2
Definition: hepunit.py:293
float proton_mass_c2
Definition: hepunit.py:274
int fine_structure_const
Definition: hepunit.py:286
float amu_c2
Definition: hepunit.py:276
int G4lrint(double ad)
Definition: templates.hh:134