21 #include "CLHEP/Units/PhysicalConstants.h" 22 #include "Randomize.hh" 23 #include "G4Electron.hh" 24 #include "G4LossTableManager.hh" 25 #include "G4ParticleChangeForMSC.hh" 27 #include "G4Poisson.hh" 33 using namespace CLHEP;
35 static const G4double
kappa = 2.5;
73 scr1ini = fine_structure_const*fine_structure_const*
74 electron_mass_c2*electron_mass_c2/(0.885*0.885*4.*
pi);
75 scr2ini = 3.76*fine_structure_const*fine_structure_const;
95 mass = proton_mass_c2;
121 if(p->GetPDGMass() >
MeV) {
122 G4cout <<
"### WARNING: UrbanMscModel93 model is used for " 123 << p->GetParticleName() <<
" !!! " << G4endl;
124 G4cout <<
"### This model should be used only for e+-" 134 const G4ParticleDefinition*
part,
135 G4double KineticEnergy,
136 G4double AtomicNumber,G4double,
139 static const G4double sigmafactor =
140 twopi*classic_electr_radius*classic_electr_radius;
141 static const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
143 static const G4double epsmin = 1.e-4 , epsmax = 1.e10;
145 static const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38., 47.,
146 50., 56., 64., 74., 79., 82. };
148 static const G4double Tdat[22] = { 100*eV, 200*eV, 400*eV, 700*eV,
149 1*keV, 2*keV, 4*keV, 7*keV,
150 10*keV, 20*keV, 40*keV, 70*keV,
151 100*keV, 200*keV, 400*keV, 700*keV,
156 static const G4double celectron[15][22] =
157 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
158 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
159 1.112,1.108,1.100,1.093,1.089,1.087 },
160 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
161 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
162 1.109,1.105,1.097,1.090,1.086,1.082 },
163 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
164 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
165 1.131,1.124,1.113,1.104,1.099,1.098 },
166 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
167 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
168 1.112,1.105,1.096,1.089,1.085,1.098 },
169 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
170 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
171 1.073,1.070,1.064,1.059,1.056,1.056 },
172 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
173 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
174 1.074,1.070,1.063,1.059,1.056,1.052 },
175 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
176 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
177 1.068,1.064,1.059,1.054,1.051,1.050 },
178 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
179 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
180 1.039,1.037,1.034,1.031,1.030,1.036 },
181 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
182 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
183 1.031,1.028,1.024,1.022,1.021,1.024 },
184 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
185 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
186 1.020,1.017,1.015,1.013,1.013,1.020 },
187 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
188 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
189 0.995,0.993,0.993,0.993,0.993,1.011 },
190 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
191 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
192 0.974,0.972,0.973,0.974,0.975,0.987 },
193 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
194 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
195 0.950,0.947,0.949,0.952,0.954,0.963 },
196 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
197 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
198 0.941,0.938,0.940,0.944,0.946,0.954 },
199 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
200 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
201 0.933,0.930,0.933,0.936,0.939,0.949 }};
203 static const G4double cpositron[15][22] = {
204 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
205 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
206 1.131,1.126,1.117,1.108,1.103,1.100 },
207 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
208 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
209 1.138,1.132,1.122,1.113,1.108,1.102 },
210 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
211 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
212 1.203,1.190,1.173,1.159,1.151,1.145 },
213 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
214 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
215 1.225,1.210,1.191,1.175,1.166,1.174 },
216 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
217 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
218 1.217,1.203,1.184,1.169,1.160,1.151 },
219 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
220 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
221 1.237,1.222,1.201,1.184,1.174,1.159 },
222 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
223 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
224 1.252,1.234,1.212,1.194,1.183,1.170 },
225 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
226 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
227 1.254,1.237,1.214,1.195,1.185,1.179 },
228 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
229 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
230 1.312,1.288,1.258,1.235,1.221,1.205 },
231 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
232 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
233 1.320,1.294,1.264,1.240,1.226,1.214 },
234 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
235 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
236 1.328,1.302,1.270,1.245,1.231,1.233 },
237 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
238 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
239 1.361,1.330,1.294,1.267,1.251,1.239 },
240 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
241 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
242 1.409,1.372,1.330,1.298,1.280,1.258 },
243 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
244 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
245 1.442,1.400,1.354,1.319,1.299,1.272 },
246 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
247 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
248 1.456,1.412,1.364,1.328,1.307,1.282 }};
251 static const G4double Tlim = 10.*
MeV;
252 static const G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
253 ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
254 static const G4double bg2lim = Tlim*(Tlim+2.*electron_mass_c2)/
255 (electron_mass_c2*electron_mass_c2);
257 static const G4double sig0[15] = {
258 0.2672*barn, 0.5922*barn, 2.653*barn, 6.235*barn,
259 11.69*barn , 13.24*barn , 16.12*barn, 23.00*barn ,
260 35.13*barn , 39.95*barn , 50.85*barn, 67.19*barn ,
261 91.15*barn , 104.4*barn , 113.1*barn};
263 static const G4double hecorr[15] = {
264 120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
265 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
271 Z23 =
pow(AtomicNumber,2./3.);
277 G4double eKineticEnergy = KineticEnergy;
279 if(
mass > electron_mass_c2)
281 G4double
TAU = KineticEnergy/
mass ;
282 G4double
c =
mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
284 G4double
tau = 0.5*(w+
sqrt(w*w+4.*c)) ;
285 eKineticEnergy = electron_mass_c2*
tau ;
288 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
289 G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
290 /(eTotalEnergy*eTotalEnergy);
291 G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
292 /(electron_mass_c2*electron_mass_c2);
294 G4double eps = epsfactor*bg2/
Z23;
296 if (eps<epsmin) sigma = 2.*eps*eps;
297 else if(eps<epsmax) sigma = G4Log(1.+2.*eps)-2.*eps/(1.+2.*eps);
298 else sigma = G4Log(2.*eps)-1.+1./eps;
300 sigma *=
ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
303 G4double
c1,c2,cc1,cc2,
corr;
307 while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
311 G4double ZZ1 = Zdat[iZ];
312 G4double ZZ2 = Zdat[iZ+1];
313 G4double ratZ = (AtomicNumber-ZZ1)*(AtomicNumber+ZZ1)/
314 ((ZZ2-ZZ1)*(ZZ2+ZZ1));
316 if(eKineticEnergy <= Tlim)
320 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
325 G4double
T = Tdat[iT], E = T + electron_mass_c2;
326 G4double b2small = T*(E+electron_mass_c2)/(E*E);
328 T = Tdat[iT+1]; E = T + electron_mass_c2;
329 G4double b2big = T*(E+electron_mass_c2)/(E*E);
330 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
334 c1 = celectron[iZ][iT];
335 c2 = celectron[iZ+1][iT];
336 cc1 = c1+ratZ*(c2-
c1);
338 c1 = celectron[iZ][iT+1];
339 c2 = celectron[iZ+1][iT+1];
340 cc2 = c1+ratZ*(c2-
c1);
342 corr = cc1+ratb2*(cc2-cc1);
344 sigma *= sigmafactor/
corr;
348 c1 = cpositron[iZ][iT];
349 c2 = cpositron[iZ+1][iT];
350 cc1 = c1+ratZ*(c2-
c1);
352 c1 = cpositron[iZ][iT+1];
353 c2 = cpositron[iZ+1][iT+1];
354 cc2 = c1+ratZ*(c2-
c1);
356 corr = cc1+ratb2*(cc2-cc1);
358 sigma *= sigmafactor/
corr;
363 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
364 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
365 if((AtomicNumber >= ZZ1) && (AtomicNumber <= ZZ2))
366 sigma = c1+ratZ*(c2-
c1) ;
367 else if(AtomicNumber < ZZ1)
368 sigma = AtomicNumber*AtomicNumber*c1/(ZZ1*ZZ1);
369 else if(AtomicNumber > ZZ2)
370 sigma = AtomicNumber*AtomicNumber*c2/(ZZ2*ZZ2);
380 SetParticle(track->GetDynamicParticle()->GetDefinition());
392 const G4Track&
track,
393 G4double& currentMinimalStep)
396 const G4DynamicParticle* dp = track.GetDynamicParticle();
397 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
398 G4StepStatus stepStatus = sp->GetStepStatus();
399 couple = track.GetMaterialCutsCouple();
408 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
424 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
429 if (steppingAlgorithm == fUseDistanceToBoundary)
438 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
444 if(
firstStep || stepStatus == fGeomBoundary)
453 rat = 1.e-3/(rat*(10.+rat)) ;
470 if(stepStatus == fGeomBoundary)
496 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
523 G4double temptlimit =
tlimit;
543 else if(steppingAlgorithm == fUseSafety)
554 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
557 if(
firstStep || stepStatus == fGeomBoundary)
572 rat = 1.e-3/(rat*(10.+rat)) ;
592 if (stepStatus == fGeomBoundary)
603 return ConvertTrueToGeom(
tPathLength, currentMinimalStep);
637 else zmean =
lambda0*(1.-G4Exp(-tau));
648 G4double lambda1 = GetTransportMeanFreePath(
particle,T1);
661 const G4double ztmax = 0.99 ;
669 G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
671 G4double u0 = cz/cz1 ;
674 u = G4Exp(G4Log(G4UniformRand())/cz1) ;
675 grej = G4Exp(cz*G4Log(u/u0))*(1.-u)/(1.-u0) ;
676 }
while (grej < G4UniformRand()) ;
681 u = 1.-G4Exp(G4Log(G4UniformRand())/cz1) ;
731 fDisplacement.set(0.0,0.0,0.0);
744 if(std::fabs(cth) > 1.) {
return fDisplacement; }
769 G4double sth =
sqrt((1.0 - cth)*(1.0 + cth));
770 G4double
phi = twopi*G4UniformRand();
771 G4double dirx = sth*
cos(phi);
772 G4double diry = sth*
sin(phi);
774 G4ThreeVector newDirection(dirx,diry,cth);
775 newDirection.rotateUz(oldDirection);
791 if(latcorr > r) latcorr =
r;
797 Phi = twopi*G4UniformRand();
800 G4double
psi = std::acos(latcorr/(r*sth));
801 if(G4UniformRand() < 0.5)
810 fDisplacement.set(dirx,diry,0.0);
811 fDisplacement.rotateUz(oldDirection);
814 return fDisplacement;
820 G4double KineticEnergy)
825 Zeff =
couple->GetMaterial()->GetTotNbOfElectPerVolume()/
826 couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
836 G4int
n = G4Poisson(mean);
840 G4double mom2 = KineticEnergy*(2.*
mass+KineticEnergy);
841 G4double beta2 = mom2/((KineticEnergy+
mass)*(KineticEnergy+
mass));
842 G4double ascr =
scr1/mom2;
843 ascr *= 1.13+
scr2/beta2;
844 G4double ascr1 = 1.+2.*ascr;
845 G4double bp1=ascr1+1.;
846 G4double bm1=ascr1-1.;
850 G4double
sx=0.,
sy=0.,sz=0.;
851 for(G4int
i=1;
i<=
n;
i++)
853 ct = ascr1-bp1*bm1/(2.*G4UniformRand()+bm1);
854 if(ct < -1.) ct = -1.;
857 phi = twopi*G4UniformRand();
869 if(
par1*trueStepLength < 1.)
870 tau = -
par2*G4Log(1.-
par1*trueStepLength) ;
880 if (tau >=
taubig) cth = -1.+2.*G4UniformRand();
883 G4double xmeanth, x2meanth;
885 xmeanth = 1.0 - tau*(1.0 - 0.5*
tau);
886 x2meanth= 1.0 - tau*(5.0 - 6.25*
tau)*
third;
888 xmeanth = G4Exp(-tau);
889 x2meanth = (1.+2.*G4Exp(-2.5*tau))*
third;
896 G4double theta0 =
ComputeTheta0(trueStepLength,KineticEnergy);
902 G4double theta2 = theta0*theta0;
904 if(theta2 <
tausmall) {
return cth; }
910 G4double
x = theta2*(1.0 - theta2/12.);
912 G4double sth = 2.*
sin(0.5*theta0);
916 G4double xmean1 = 1.-(1.-(1.+
xsi)*
ea)*x/
eaa;
917 G4double x0 = 1. -
xsi*
x;
921 if(xmean1 <= 0.999*xmeanth) {
942 if(fabs(c-3.) < 0.001) { c = 3.001; }
943 else if(fabs(c-2.) < 0.001) { c = 2.001; }
948 G4double
b = 1.+(c-
xsi)*x;
953 G4double eb1 =
pow(b1,c1);
954 G4double ebx =
pow(bx,c1);
955 G4double
d = ebx/eb1;
958 G4double xmean2 = (x0 + d - (bx - b1*
d)/(c-2.))/(1. - d);
960 G4double f1x0 =
ea/
eaa;
961 G4double f2x0 = c1/(c*(1. -
d));
962 G4double
prob = f2x0/(f1x0+f2x0);
964 G4double qprob = xmeanth/(prob*xmean1+(1.-
prob)*xmean2);
970 if(G4UniformRand() < qprob)
973 if(G4UniformRand() <
prob) {
974 cth = 1.+G4Log(
ea+G4UniformRand()*
eaa)*
x;
976 var = (1.0 -
d)*G4UniformRand();
979 cth = -1.0 + var*(1.0 - 0.5*var*
c)*(2. + (c -
xsi)*
x);
981 cth = 1. + x*(c -
xsi - c*
pow(var + d, -1.0/c1));
985 if(KineticEnergy > 5*
GeV && cth < 0.9) {
986 G4cout <<
"UrbanMscModel93::SampleCosineTheta: E(GeV)= " 988 <<
" 1-cosT= " << 1 - cth
989 <<
" length(mm)= " << trueStepLength <<
" Zeff= " <<
Zeff 991 <<
" prob= " << prob <<
" var= " << var << G4endl;
992 G4cout <<
" c= " << c <<
" qprob= " << qprob <<
" eb1= " << eb1
994 <<
" c1= " << c1 <<
" b= " << b <<
" b1= " << b1
995 <<
" bx= " << bx <<
" d= " << d
996 <<
" ea= " <<
ea <<
" eaa= " <<
eaa << G4endl;
1000 cth = -1.+2.*G4UniformRand();
1001 if(KineticEnergy > 5*
GeV) {
1002 G4cout <<
"UrbanMscModel93::SampleCosineTheta: E(GeV)= " 1003 << KineticEnergy/
GeV 1004 <<
" length(mm)= " << trueStepLength <<
" Zeff= " <<
Zeff 1005 <<
" qprob= " << qprob << G4endl;
1018 G4double rmean = 0.0;
1025 G4double etau = 0.0;
1035 if(rmean == 0.)
return rmean;
1044 if(rmean >= rmax)
return rmax;
1067 G4double latcorr = 0.;
const G4ParticleDefinition * particle
void StartTracking(G4Track *) override
const G4MaterialCutsCouple * couple
G4ParticleChangeForMSC * fParticleChange
G4int currentMaterialIndex
static const G4double kappami1
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
Sin< T >::type sin(const T &t)
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double currentRadLength
void SetParticle(const G4ParticleDefinition *)
std::map< std::string, int, std::less< std::string > > psi
G4double ComputeGeomPathLength(G4double truePathLength) override
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double ¤tMinimalStep) override
G4double LatCorrelation()
UrbanMscModel93(const G4String &nam="UrbanMsc93")
G4LossTableManager * theManager
G4double currentKinEnergy
Cos< T >::type cos(const T &t)
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
G4double SampleDisplacement()
Abs< T >::type abs(const T &t)
~UrbanMscModel93() override
static const G4double kappapl1
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeTrueStepLength(G4double geomStepLength) override
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
static const G4double kappa
Power< A, B >::type pow(const A &a, const B &b)