51 #include "CMSUrbanMscModel92.hh"
52 #include "Randomize.hh"
53 #include "G4Electron.hh"
54 #include "G4LossTableManager.hh"
55 #include "G4ParticleChangeForMSC.hh"
57 #include "G4Poisson.hh"
64 CMSUrbanMscModel92::CMSUrbanMscModel92(
const G4String& nam)
76 tlimitminfix = 1.e-6*mm;
77 stepmin = tlimitminfix;
82 tlimitmin = 10.*tlimitminfix;
100 scr1ini = fine_structure_const*fine_structure_const*
101 electron_mass_c2*electron_mass_c2/(0.885*0.885*4.*
pi);
102 scr2ini = 3.76*fine_structure_const*fine_structure_const;
110 theManager = G4LossTableManager::Instance();
114 skindepth = skin*stepmin;
116 mass = proton_mass_c2;
117 charge = ChargeSquare = 1.0;
118 currentKinEnergy = currentRadLength = lambda0 = lambdaeff = tPathLength
119 = zPathLength = par1 = par2 = par3 = 0;
121 currentMaterialIndex = 0;
126 CMSUrbanMscModel92::~CMSUrbanMscModel92()
131 void CMSUrbanMscModel92::Initialise(
const G4ParticleDefinition*
p,
134 skindepth = skin*stepmin;
135 if(isInitialized)
return;
139 if(p->GetPDGMass() > MeV) {
140 G4cout <<
"### WARNING: CMSUrbanMscModel92 model is used for "
141 << p->GetParticleName() <<
" !!! " << G4endl;
142 G4cout <<
"### This model should be used only for e+-"
146 fParticleChange = GetParticleChangeForMSC();
147 InitialiseSafetyHelper();
149 isInitialized =
true;
154 G4double CMSUrbanMscModel92::ComputeCrossSectionPerAtom(
155 const G4ParticleDefinition*
part,
156 G4double KineticEnergy,
157 G4double AtomicNumber,G4double,
160 const G4double sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
161 const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
162 Bohr_radius*Bohr_radius/(hbarc*hbarc);
163 const G4double epsmin = 1.e-4 , epsmax = 1.e10;
165 const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38., 47.,
166 50., 56., 64., 74., 79., 82. };
168 const G4double Tdat[22] = { 100*eV, 200*eV, 400*eV, 700*eV,
169 1*keV, 2*keV, 4*keV, 7*keV,
170 10*keV, 20*keV, 40*keV, 70*keV,
171 100*keV, 200*keV, 400*keV, 700*keV,
172 1*MeV, 2*MeV, 4*MeV, 7*MeV,
176 G4double celectron[15][22] =
177 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
178 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
179 1.112,1.108,1.100,1.093,1.089,1.087 },
180 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
181 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
182 1.109,1.105,1.097,1.090,1.086,1.082 },
183 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
184 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
185 1.131,1.124,1.113,1.104,1.099,1.098 },
186 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
187 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
188 1.112,1.105,1.096,1.089,1.085,1.098 },
189 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
190 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
191 1.073,1.070,1.064,1.059,1.056,1.056 },
192 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
193 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
194 1.074,1.070,1.063,1.059,1.056,1.052 },
195 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
196 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
197 1.068,1.064,1.059,1.054,1.051,1.050 },
198 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
199 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
200 1.039,1.037,1.034,1.031,1.030,1.036 },
201 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
202 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
203 1.031,1.028,1.024,1.022,1.021,1.024 },
204 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
205 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
206 1.020,1.017,1.015,1.013,1.013,1.020 },
207 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
208 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
209 0.995,0.993,0.993,0.993,0.993,1.011 },
210 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
211 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
212 0.974,0.972,0.973,0.974,0.975,0.987 },
213 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
214 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
215 0.950,0.947,0.949,0.952,0.954,0.963 },
216 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
217 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
218 0.941,0.938,0.940,0.944,0.946,0.954 },
219 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
220 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
221 0.933,0.930,0.933,0.936,0.939,0.949 }};
223 G4double cpositron[15][22] = {
224 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
225 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
226 1.131,1.126,1.117,1.108,1.103,1.100 },
227 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
228 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
229 1.138,1.132,1.122,1.113,1.108,1.102 },
230 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
231 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
232 1.203,1.190,1.173,1.159,1.151,1.145 },
233 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
234 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
235 1.225,1.210,1.191,1.175,1.166,1.174 },
236 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
237 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
238 1.217,1.203,1.184,1.169,1.160,1.151 },
239 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
240 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
241 1.237,1.222,1.201,1.184,1.174,1.159 },
242 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
243 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
244 1.252,1.234,1.212,1.194,1.183,1.170 },
245 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
246 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
247 1.254,1.237,1.214,1.195,1.185,1.179 },
248 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
249 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
250 1.312,1.288,1.258,1.235,1.221,1.205 },
251 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
252 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
253 1.320,1.294,1.264,1.240,1.226,1.214 },
254 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
255 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
256 1.328,1.302,1.270,1.245,1.231,1.233 },
257 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
258 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
259 1.361,1.330,1.294,1.267,1.251,1.239 },
260 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
261 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
262 1.409,1.372,1.330,1.298,1.280,1.258 },
263 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
264 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
265 1.442,1.400,1.354,1.319,1.299,1.272 },
266 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
267 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
268 1.456,1.412,1.364,1.328,1.307,1.282 }};
271 G4double Tlim = 10.*MeV;
272 G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
273 ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
274 G4double bg2lim = Tlim*(Tlim+2.*electron_mass_c2)/
275 (electron_mass_c2*electron_mass_c2);
277 G4double sig0[15] = {0.2672*barn, 0.5922*barn, 2.653*barn, 6.235*barn,
278 11.69*barn , 13.24*barn , 16.12*barn, 23.00*barn ,
279 35.13*barn , 39.95*barn , 50.85*barn, 67.19*barn ,
280 91.15*barn , 104.4*barn , 113.1*barn};
282 G4double hecorr[15] = {120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
283 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
289 G4double Z23 = 2.*
log(AtomicNumber)/3.; Z23 =
exp(Z23);
295 G4double eKineticEnergy = KineticEnergy;
297 if(
mass > electron_mass_c2)
299 G4double
TAU = KineticEnergy/
mass ;
300 G4double
c =
mass*TAU*(TAU+2.)/(electron_mass_c2*(
TAU+1.)) ;
302 G4double
tau = 0.5*(w+
sqrt(w*w+4.*
c)) ;
303 eKineticEnergy = electron_mass_c2*
tau ;
306 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
307 G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
308 /(eTotalEnergy*eTotalEnergy);
309 G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
310 /(electron_mass_c2*electron_mass_c2);
312 G4double eps = epsfactor*bg2/Z23;
314 if (eps<epsmin) sigma = 2.*eps*eps;
315 else if(eps<epsmax) sigma =
log(1.+2.*eps)-2.*eps/(1.+2.*eps);
316 else sigma =
log(2.*eps)-1.+1./eps;
318 sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
321 G4double
c1,c2,cc1,cc2,
corr;
325 while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
329 G4double Z1 = Zdat[iZ];
330 G4double Z2 = Zdat[iZ+1];
331 G4double ratZ = (AtomicNumber-Z1)*(AtomicNumber+Z1)/
334 if(eKineticEnergy <= Tlim)
338 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
343 G4double
T = Tdat[iT], E = T + electron_mass_c2;
344 G4double b2small = T*(E+electron_mass_c2)/(E*E);
346 T = Tdat[iT+1]; E = T + electron_mass_c2;
347 G4double b2big = T*(E+electron_mass_c2)/(E*E);
348 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
352 c1 = celectron[iZ][iT];
353 c2 = celectron[iZ+1][iT];
354 cc1 = c1+ratZ*(c2-
c1);
356 c1 = celectron[iZ][iT+1];
357 c2 = celectron[iZ+1][iT+1];
358 cc2 = c1+ratZ*(c2-
c1);
360 corr = cc1+ratb2*(cc2-cc1);
362 sigma *= sigmafactor/
corr;
366 c1 = cpositron[iZ][iT];
367 c2 = cpositron[iZ+1][iT];
368 cc1 = c1+ratZ*(c2-
c1);
370 c1 = cpositron[iZ][iT+1];
371 c2 = cpositron[iZ+1][iT+1];
372 cc2 = c1+ratZ*(c2-
c1);
374 corr = cc1+ratb2*(cc2-cc1);
376 sigma *= sigmafactor/
corr;
381 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
382 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
383 if((AtomicNumber >= Z1) && (AtomicNumber <= Z2))
384 sigma = c1+ratZ*(c2-
c1) ;
385 else if(AtomicNumber < Z1)
386 sigma = AtomicNumber*AtomicNumber*c1/(Z1*Z1);
387 else if(AtomicNumber > Z2)
388 sigma = AtomicNumber*AtomicNumber*c2/(Z2*Z2);
396 G4double CMSUrbanMscModel92::ComputeTruePathLengthLimit(
397 const G4Track& track,
398 G4PhysicsTable* theTable,
399 G4double currentMinimalStep)
401 tPathLength = currentMinimalStep;
402 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
403 G4StepStatus stepStatus = sp->GetStepStatus();
405 const G4DynamicParticle* dp = track.GetDynamicParticle();
407 if(stepStatus == fUndefined) {
411 SetParticle( dp->GetDefinition() );
414 theLambdaTable = theTable;
415 couple = track.GetMaterialCutsCouple();
416 currentMaterialIndex = couple->GetIndex();
417 currentKinEnergy = dp->GetKineticEnergy();
419 theManager->GetRangeFromRestricteDEDX(particle,currentKinEnergy,couple);
420 lambda0 = GetLambda(currentKinEnergy);
423 if(inside)
return tPathLength;
425 if(tPathLength > currentRange) tPathLength = currentRange;
427 presafety = sp->GetSafety();
434 if(currentRange < presafety)
442 if (steppingAlgorithm == fUseDistanceToBoundary)
445 G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
448 if(currentRange < presafety)
457 if((stepStatus == fGeomBoundary) || (stepStatus == fUndefined))
459 rangeinit = currentRange;
460 if(stepStatus == fUndefined) smallstep = 1.e10;
464 if((geomlimit < geombig) && (geomlimit > geommin))
466 if(stepStatus == fGeomBoundary)
467 tgeom = geomlimit/facgeom;
469 tgeom = 2.*geomlimit/facgeom;
476 G4double rat = currentKinEnergy/MeV ;
477 rat = 1.e-3/(rat*(10.+rat)) ;
479 stepmin = rat*lambda0;
480 skindepth = skin*stepmin;
483 tlimitmin = 10.*stepmin;
484 if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
489 tlimit = facrange*rangeinit;
490 if(tlimit < facsafety*presafety)
491 tlimit = facsafety*presafety;
494 if(tlimit < tlimitmin) tlimit = tlimitmin;
496 if(tlimit > tgeom) tlimit = tgeom;
502 if((tPathLength < tlimit) && (tPathLength < presafety) &&
503 (smallstep >= skin) && (tPathLength < geomlimit-0.999*skindepth))
512 else if(geomlimit < geombig)
514 if(geomlimit > skindepth)
516 if(tlimit > geomlimit-0.999*skindepth)
517 tlimit = geomlimit-0.999*skindepth;
522 if(tlimit > stepmin) tlimit = stepmin;
526 if(tlimit < stepmin) tlimit = stepmin;
528 if(tPathLength > tlimit) tPathLength = tlimit ;
533 else if(steppingAlgorithm == fUseSafety)
537 if((stepStatus != fGeomBoundary) && (presafety < tlimitminfix))
538 presafety = ComputeSafety(sp->GetPosition(),tPathLength);
541 if(currentRange < presafety)
547 if((stepStatus == fGeomBoundary) || (stepStatus == fUndefined))
549 rangeinit = currentRange;
552 if(
mass < masslimite)
554 if(lambda0 > currentRange)
556 if(lambda0 > lambdalimit)
557 fr *= 0.75+0.25*lambda0/lambdalimit;
561 G4double rat = currentKinEnergy/MeV ;
562 rat = 1.e-3/(rat*(10.+rat)) ;
563 tlimitmin = 10.*lambda0*rat;
564 if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
567 tlimit = fr*rangeinit;
569 if(tlimit < facsafety*presafety)
570 tlimit = facsafety*presafety;
573 if(tlimit < tlimitmin) tlimit = tlimitmin;
575 if(tPathLength > tlimit) tPathLength = tlimit;
581 if (stepStatus == fGeomBoundary)
583 if (currentRange > lambda0) tlimit = facrange*currentRange;
584 else tlimit = facrange*lambda0;
586 if(tlimit < tlimitmin) tlimit = tlimitmin;
587 if(tPathLength > tlimit) tPathLength = tlimit;
597 G4double CMSUrbanMscModel92::ComputeGeomPathLength(G4double)
604 zPathLength = tPathLength;
607 if(tPathLength < tlimitminfix)
return zPathLength;
611 if(tPathLength > currentRange)
612 tPathLength = currentRange ;
614 G4double
tau = tPathLength/lambda0 ;
616 if ((tau <= tausmall) || insideskin) {
617 zPathLength = tPathLength;
618 if(zPathLength > lambda0) zPathLength = lambda0;
622 G4double zmean = tPathLength;
623 if (tPathLength < currentRange*dtrl) {
624 if(tau < taulim) zmean = tPathLength*(1.-0.5*
tau) ;
625 else zmean = lambda0*(1.-
exp(-tau));
626 }
else if(currentKinEnergy <
mass ) {
627 par1 = 1./currentRange ;
628 par2 = 1./(par1*lambda0) ;
630 if(tPathLength < currentRange)
631 zmean = (1.-
exp(par3*
log(1.-tPathLength/currentRange)))/(par1*par3) ;
633 zmean = 1./(par1*par3) ;
635 G4double T1 = theManager->GetEnergy(particle,currentRange-tPathLength,couple);
636 G4double lambda1 = GetLambda(T1);
638 par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
639 par2 = 1./(par1*lambda0) ;
641 zmean = (1.-
exp(par3*
log(lambda1/lambda0)))/(par1*par3) ;
644 zPathLength = zmean ;
649 const G4double ztmax = 0.99 ;
650 G4double zt = zmean/tPathLength ;
652 if (tPathLength > stepmin && zt < ztmax)
657 G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
659 G4double u0 = cz/cz1 ;
662 u =
exp(
log(G4UniformRand())/cz1) ;
663 grej =
exp(cz*
log(u/u0))*(1.-u)/(1.-u0) ;
664 }
while (grej < G4UniformRand()) ;
669 u = 1.-
exp(
log(G4UniformRand())/cz1) ;
671 zPathLength = tPathLength*u ;
675 if(zPathLength > lambda0) zPathLength = lambda0;
682 G4double CMSUrbanMscModel92::ComputeTrueStepLength(G4double geomStepLength)
685 if(geomStepLength == zPathLength && tPathLength <= currentRange)
689 zPathLength = geomStepLength;
690 tPathLength = geomStepLength;
691 if(geomStepLength < tlimitminfix)
return tPathLength;
694 if((geomStepLength > lambda0*tausmall) && !insideskin)
697 tPathLength = -lambda0*
log(1.-geomStepLength/lambda0) ;
700 if(par1*par3*geomStepLength < 1.)
701 tPathLength = (1.-
exp(
log(1.-par1*par3*geomStepLength)/par3))/par1 ;
703 tPathLength = currentRange;
706 if(tPathLength < geomStepLength) tPathLength = geomStepLength;
713 G4double CMSUrbanMscModel92::ComputeTheta0(G4double trueStepLength,
714 G4double KineticEnergy)
719 const G4double c_highland = 13.6*MeV ;
720 G4double betacp =
sqrt(currentKinEnergy*(currentKinEnergy+2.*
mass)*
721 KineticEnergy*(KineticEnergy+2.*
mass)/
722 ((currentKinEnergy+
mass)*(KineticEnergy+
mass)));
723 y = trueStepLength/currentRadLength;
727 G4double corr = coeffth1+coeffth2*
y;
728 if(y < -6.5) corr -= 0.011*(6.5+
y);
736 void CMSUrbanMscModel92::SampleScattering(
const G4DynamicParticle* dynParticle,
739 G4double kineticEnergy = dynParticle->GetKineticEnergy();
741 if((kineticEnergy <= 0.0) || (tPathLength <= tlimitminfix) ||
742 (tPathLength/tausmall < lambda0))
return;
744 G4double cth = SampleCosineTheta(tPathLength,kineticEnergy);
750 if(cth < 1.0 - 1000*tPathLength/lambda0 && kineticEnergy > 20*MeV) {
756 cth = 1.0 + 2*
log(G4UniformRand())*tPathLength/lambda0;
760 G4double sth =
sqrt((1.0 - cth)*(1.0 + cth));
761 G4double
phi = twopi*G4UniformRand();
762 G4double dirx = sth*
cos(phi);
763 G4double diry = sth*
sin(phi);
765 G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
766 G4ThreeVector newDirection(dirx,diry,cth);
767 newDirection.rotateUz(oldDirection);
768 fParticleChange->ProposeMomentumDirection(newDirection);
770 if (latDisplasment && safety > tlimitminfix) {
772 G4double
r = SampleDisplacement();
782 G4double latcorr = LatCorrelation();
783 if(latcorr > r) latcorr =
r;
789 Phi = twopi*G4UniformRand();
792 G4double
psi = std::acos(latcorr/(r*sth));
793 if(G4UniformRand() < 0.5)
802 G4ThreeVector latDirection(dirx,diry,0.0);
803 latDirection.rotateUz(oldDirection);
805 ComputeDisplacement(fParticleChange, latDirection, r, safety);
812 G4double CMSUrbanMscModel92::SampleCosineTheta(G4double trueStepLength,
813 G4double KineticEnergy)
816 G4double tau = trueStepLength/lambda0 ;
818 Zeff = couple->GetMaterial()->GetTotNbOfElectPerVolume()/
819 couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
827 G4double
mean = trueStepLength/stepmin ;
829 G4int
n = G4Poisson(mean);
833 G4double mom2 = KineticEnergy*(2.*
mass+KineticEnergy);
834 G4double beta2 = mom2/((KineticEnergy+
mass)*(KineticEnergy+
mass));
835 G4double ascr = scr1/mom2;
836 ascr *= 1.13+scr2/beta2;
837 G4double ascr1 = 1.+2.*ascr;
838 G4double bp1=ascr1+1.;
839 G4double bm1=ascr1-1.;
843 G4double sx=0.,sy=0.,sz=0.;
844 for(G4int
i=1;
i<=
n;
i++)
846 ct = ascr1-bp1*bm1/(2.*G4UniformRand()+bm1);
847 if(ct < -1.) ct = -1.;
850 phi = twopi*G4UniformRand();
855 cth = sz/
sqrt(sx*sx+sy*sy+sz*sz);
860 if(trueStepLength >= currentRange*dtrl)
862 if(par1*trueStepLength < 1.)
863 tau = -par2*
log(1.-par1*trueStepLength) ;
866 else if(1.-KineticEnergy/currentKinEnergy > taulim)
870 lambdaeff = trueStepLength/currentTau;
871 currentRadLength = couple->GetMaterial()->GetRadlen();
873 if (tau >= taubig) cth = -1.+2.*G4UniformRand();
874 else if (tau >= tausmall)
878 G4double
a = 1., ea = 0., eaa = 1.;
879 G4double
b=2.,b1=3.,bx=1.,eb1=3.,ebx=1.;
880 G4double prob = 1. , qprob = 1. ;
881 G4double xmean1 = 1., xmean2 = 0.;
882 G4double xmeanth =
exp(-tau);
883 G4double x2meanth = (1.+2.*
exp(-2.5*tau))/3.;
885 G4double relloss = 1.-KineticEnergy/currentKinEnergy;
886 if(relloss > rellossmax)
887 return SimpleScattering(xmeanth,x2meanth);
889 G4double theta0 = ComputeTheta0(trueStepLength,KineticEnergy);
892 if(theta0*theta0 < tausmall)
return cth;
894 if(theta0 > theta0max)
895 return SimpleScattering(xmeanth,x2meanth);
896 G4double sth =
sin(0.5*theta0);
901 xmean1 = 1.-(1.-(1.+xsi)*ea)/(a*eaa);
904 if(xmean1 <= 0.999*xmeanth)
905 return SimpleScattering(xmeanth,x2meanth);
908 G4double
c = coeffc1;
910 c += coeffc2*
exp(3.*
log(y+13.5));
912 if(
abs(c-3.) < 0.001) c = 3.001;
913 if(
abs(c-2.) < 0.001) c = 2.001;
925 xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/(c-2.))/(eb1-ebx);
927 G4double f1x0 = a*ea/eaa;
928 G4double f2x0 = c1*eb1/(bx*(eb1-ebx));
929 prob = f2x0/(f1x0+f2x0);
931 qprob = xmeanth/(prob*xmean1+(1.-prob)*xmean2);
934 if(G4UniformRand() < qprob)
936 if(G4UniformRand() < prob)
937 cth = 1.+
log(ea+G4UniformRand()*eaa)/a ;
939 cth = b-b1*bx/
exp(
log(ebx+(eb1-ebx)*G4UniformRand())/c1) ;
942 cth = -1.+2.*G4UniformRand();
950 G4double CMSUrbanMscModel92::SimpleScattering(G4double xmeanth,G4double x2meanth)
954 G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
955 G4double prob = (a+2.)*xmeanth/
a;
959 if(G4UniformRand() < prob)
960 cth = -1.+2.*
exp(
log(G4UniformRand())/(a+1.));
962 cth = -1.+2.*G4UniformRand();
968 G4double CMSUrbanMscModel92::SampleDisplacement()
970 const G4double kappa = 2.5;
971 const G4double kappapl1 = kappa+1.;
972 const G4double kappami1 = kappa-1.;
973 G4double rmean = 0.0;
974 if ((currentTau >= tausmall) && !insideskin) {
975 if (currentTau < taulim) {
976 rmean = kappa*currentTau*currentTau*currentTau*
977 (1.-kappapl1*currentTau*0.25)/6. ;
981 if (currentTau<taubig) etau =
exp(-currentTau);
982 rmean = -kappa*currentTau;
983 rmean = -
exp(rmean)/(kappa*kappami1);
984 rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
986 if (rmean>0.) rmean = 2.*lambdaeff*
sqrt(rmean/3.0);
992 G4double zt = (tPathLength-zPathLength)*(tPathLength+zPathLength);
995 else if(rmean*rmean > zt)
1003 G4double CMSUrbanMscModel92::LatCorrelation()
1005 const G4double kappa = 2.5;
1006 const G4double kappami1 = kappa-1.;
1008 G4double latcorr = 0.;
1009 if((currentTau >= tausmall) && !insideskin)
1011 if(currentTau < taulim)
1012 latcorr = lambdaeff*kappa*currentTau*currentTau*
1013 (1.-(kappa+1.)*currentTau/3.)/3.;
1017 if(currentTau < taubig) etau =
exp(-currentTau);
1018 latcorr = -kappa*currentTau;
1019 latcorr =
exp(latcorr)/kappami1;
1020 latcorr += 1.-kappa*etau/kappami1 ;
1021 latcorr *= 2.*lambdaeff/3. ;
Sin< T >::type sin(const T &t)
std::map< std::string, int, std::less< std::string > > psi
Cos< T >::type cos(const T &t)