00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056 #include "G4UrbanMscModel90.hh"
00057 #include "Randomize.hh"
00058 #include "G4Electron.hh"
00059 #include "G4LossTableManager.hh"
00060 #include "G4ParticleChangeForMSC.hh"
00061 #include "G4TransportationManager.hh"
00062 #include "G4SafetyHelper.hh"
00063
00064 #include "G4Poisson.hh"
00065
00066
00067
00068 using namespace std;
00069
00070 G4UrbanMscModel90::G4UrbanMscModel90(G4double m_facrange, G4double m_dtrl,
00071 G4double m_lambdalimit,
00072 G4double m_facgeom,G4double m_skin,
00073 G4bool m_samplez, G4MscStepLimitType m_stepAlg,
00074 const G4String& nam)
00075 : G4VEmModel(nam),
00076 dtrl(m_dtrl),
00077 lambdalimit(m_lambdalimit),
00078 facrange(m_facrange),
00079 facgeom(m_facgeom),
00080 skin(m_skin),
00081 steppingAlgorithm(m_stepAlg),
00082 samplez(m_samplez),
00083 isInitialized(false)
00084 {
00085 taubig = 8.0;
00086 tausmall = 1.e-20;
00087 taulim = 1.e-6;
00088 currentTau = taulim;
00089 tlimitminfix = 1.e-6*mm;
00090 stepmin = tlimitminfix;
00091 skindepth = skin*stepmin;
00092 smallstep = 1.e10;
00093 currentRange = 0. ;
00094 frscaling2 = 0.25;
00095 frscaling1 = 1.-frscaling2;
00096 tlimit = 1.e10*mm;
00097 tlimitmin = 10.*tlimitminfix;
00098 nstepmax = 25.;
00099 geombig = 1.e50*mm;
00100 geommin = 1.e-3*mm;
00101 geomlimit = geombig;
00102 presafety = 0.*mm;
00103 facsafety = 0.25;
00104 Zeff = 1.;
00105 particle = 0;
00106 theManager = G4LossTableManager::Instance();
00107 inside = false;
00108 insideskin = false;
00109
00110 }
00111
00112
00113
00114 G4UrbanMscModel90::~G4UrbanMscModel90()
00115 {}
00116
00117
00118
00119 void G4UrbanMscModel90::Initialise(const G4ParticleDefinition* p,
00120 const G4DataVector&)
00121 {
00122 if(isInitialized) return;
00123
00124 SetParticle(p);
00125
00126 if (pParticleChange)
00127 fParticleChange = reinterpret_cast<G4ParticleChangeForMSC*>(pParticleChange);
00128 else
00129 fParticleChange = new G4ParticleChangeForMSC();
00130
00131 safetyHelper = G4TransportationManager::GetTransportationManager()
00132 ->GetSafetyHelper();
00133 safetyHelper->InitialiseHelper();
00134 }
00135
00136
00137
00138 G4double G4UrbanMscModel90::ComputeCrossSectionPerAtom(
00139 const G4ParticleDefinition* part,
00140 G4double KineticEnergy,
00141 G4double AtomicNumber,G4double,
00142 G4double, G4double)
00143 {
00144 const G4double sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
00145 const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
00146 Bohr_radius*Bohr_radius/(hbarc*hbarc);
00147 const G4double epsmin = 1.e-4 , epsmax = 1.e10;
00148
00149 const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38., 47.,
00150 50., 56., 64., 74., 79., 82. };
00151
00152 const G4double Tdat[22] = { 100*eV, 200*eV, 400*eV, 700*eV,
00153 1*keV, 2*keV, 4*keV, 7*keV,
00154 10*keV, 20*keV, 40*keV, 70*keV,
00155 100*keV, 200*keV, 400*keV, 700*keV,
00156 1*MeV, 2*MeV, 4*MeV, 7*MeV,
00157 10*MeV, 20*MeV};
00158
00159
00160 G4double celectron[15][22] =
00161 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
00162 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
00163 1.112,1.108,1.100,1.093,1.089,1.087 },
00164 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
00165 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
00166 1.109,1.105,1.097,1.090,1.086,1.082 },
00167 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
00168 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
00169 1.131,1.124,1.113,1.104,1.099,1.098 },
00170 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
00171 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
00172 1.112,1.105,1.096,1.089,1.085,1.098 },
00173 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
00174 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
00175 1.073,1.070,1.064,1.059,1.056,1.056 },
00176 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
00177 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
00178 1.074,1.070,1.063,1.059,1.056,1.052 },
00179 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
00180 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
00181 1.068,1.064,1.059,1.054,1.051,1.050 },
00182 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
00183 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
00184 1.039,1.037,1.034,1.031,1.030,1.036 },
00185 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
00186 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
00187 1.031,1.028,1.024,1.022,1.021,1.024 },
00188 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
00189 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
00190 1.020,1.017,1.015,1.013,1.013,1.020 },
00191 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
00192 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
00193 0.995,0.993,0.993,0.993,0.993,1.011 },
00194 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
00195 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
00196 0.974,0.972,0.973,0.974,0.975,0.987 },
00197 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
00198 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
00199 0.950,0.947,0.949,0.952,0.954,0.963 },
00200 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
00201 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
00202 0.941,0.938,0.940,0.944,0.946,0.954 },
00203 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
00204 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
00205 0.933,0.930,0.933,0.936,0.939,0.949 }};
00206
00207 G4double cpositron[15][22] = {
00208 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
00209 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
00210 1.131,1.126,1.117,1.108,1.103,1.100 },
00211 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
00212 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
00213 1.138,1.132,1.122,1.113,1.108,1.102 },
00214 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
00215 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
00216 1.203,1.190,1.173,1.159,1.151,1.145 },
00217 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
00218 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
00219 1.225,1.210,1.191,1.175,1.166,1.174 },
00220 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
00221 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
00222 1.217,1.203,1.184,1.169,1.160,1.151 },
00223 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
00224 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
00225 1.237,1.222,1.201,1.184,1.174,1.159 },
00226 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
00227 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
00228 1.252,1.234,1.212,1.194,1.183,1.170 },
00229 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
00230 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
00231 1.254,1.237,1.214,1.195,1.185,1.179 },
00232 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
00233 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
00234 1.312,1.288,1.258,1.235,1.221,1.205 },
00235 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
00236 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
00237 1.320,1.294,1.264,1.240,1.226,1.214 },
00238 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
00239 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
00240 1.328,1.302,1.270,1.245,1.231,1.233 },
00241 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
00242 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
00243 1.361,1.330,1.294,1.267,1.251,1.239 },
00244 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
00245 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
00246 1.409,1.372,1.330,1.298,1.280,1.258 },
00247 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
00248 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
00249 1.442,1.400,1.354,1.319,1.299,1.272 },
00250 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
00251 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
00252 1.456,1.412,1.364,1.328,1.307,1.282 }};
00253
00254
00255 G4double Tlim = 10.*MeV;
00256 G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
00257 ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
00258 G4double bg2lim = Tlim*(Tlim+2.*electron_mass_c2)/
00259 (electron_mass_c2*electron_mass_c2);
00260
00261 G4double sig0[15] = {0.2672*barn, 0.5922*barn, 2.653*barn, 6.235*barn,
00262 11.69*barn , 13.24*barn , 16.12*barn, 23.00*barn ,
00263 35.13*barn , 39.95*barn , 50.85*barn, 67.19*barn ,
00264 91.15*barn , 104.4*barn , 113.1*barn};
00265
00266 G4double hecorr[15] = {120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
00267 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
00268 -22.30};
00269
00270 G4double sigma;
00271 SetParticle(part);
00272
00273 G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
00274
00275
00276
00277
00278
00279 G4double eKineticEnergy = KineticEnergy;
00280
00281 if((particle->GetParticleName() != "e-") &&
00282 (particle->GetParticleName() != "e+") )
00283 {
00284 G4double TAU = KineticEnergy/mass ;
00285 G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
00286 G4double w = c-2. ;
00287 G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
00288 eKineticEnergy = electron_mass_c2*tau ;
00289 }
00290
00291 G4double ChargeSquare = charge*charge;
00292
00293 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
00294 G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00295 /(eTotalEnergy*eTotalEnergy);
00296 G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00297 /(electron_mass_c2*electron_mass_c2);
00298
00299 G4double eps = epsfactor*bg2/Z23;
00300
00301 if (eps<epsmin) sigma = 2.*eps*eps;
00302 else if(eps<epsmax) sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
00303 else sigma = log(2.*eps)-1.+1./eps;
00304
00305 sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
00306
00307
00308 G4double c1,c2,cc1,cc2,corr;
00309
00310
00311 G4int iZ = 14;
00312 while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
00313 if (iZ==14) iZ = 13;
00314 if (iZ==-1) iZ = 0 ;
00315
00316 G4double Z1 = Zdat[iZ];
00317 G4double Z2 = Zdat[iZ+1];
00318 G4double ratZ = (AtomicNumber-Z1)*(AtomicNumber+Z1)/
00319 ((Z2-Z1)*(Z2+Z1));
00320
00321 if(eKineticEnergy <= Tlim)
00322 {
00323
00324 G4int iT = 21;
00325 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
00326 if(iT==21) iT = 20;
00327 if(iT==-1) iT = 0 ;
00328
00329
00330 G4double T = Tdat[iT], E = T + electron_mass_c2;
00331 G4double b2small = T*(E+electron_mass_c2)/(E*E);
00332
00333 T = Tdat[iT+1]; E = T + electron_mass_c2;
00334 G4double b2big = T*(E+electron_mass_c2)/(E*E);
00335 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
00336
00337 if (charge < 0.)
00338 {
00339 c1 = celectron[iZ][iT];
00340 c2 = celectron[iZ+1][iT];
00341 cc1 = c1+ratZ*(c2-c1);
00342
00343 c1 = celectron[iZ][iT+1];
00344 c2 = celectron[iZ+1][iT+1];
00345 cc2 = c1+ratZ*(c2-c1);
00346
00347 corr = cc1+ratb2*(cc2-cc1);
00348
00349 sigma *= sigmafactor/corr;
00350 }
00351 else
00352 {
00353 c1 = cpositron[iZ][iT];
00354 c2 = cpositron[iZ+1][iT];
00355 cc1 = c1+ratZ*(c2-c1);
00356
00357 c1 = cpositron[iZ][iT+1];
00358 c2 = cpositron[iZ+1][iT+1];
00359 cc2 = c1+ratZ*(c2-c1);
00360
00361 corr = cc1+ratb2*(cc2-cc1);
00362
00363 sigma *= sigmafactor/corr;
00364 }
00365 }
00366 else
00367 {
00368 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
00369 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
00370 if((AtomicNumber >= Z1) && (AtomicNumber <= Z2))
00371 sigma = c1+ratZ*(c2-c1) ;
00372 else if(AtomicNumber < Z1)
00373 sigma = AtomicNumber*AtomicNumber*c1/(Z1*Z1);
00374 else if(AtomicNumber > Z2)
00375 sigma = AtomicNumber*AtomicNumber*c2/(Z2*Z2);
00376 }
00377 return sigma;
00378
00379 }
00380
00381
00382
00383 G4double G4UrbanMscModel90::ComputeTruePathLengthLimit(
00384 const G4Track& track,
00385 G4PhysicsTable* theTable,
00386 G4double currentMinimalStep)
00387 {
00388 tPathLength = currentMinimalStep;
00389 G4int stepNumber = track.GetCurrentStepNumber();
00390 const G4DynamicParticle* dp = track.GetDynamicParticle();
00391
00392 if(stepNumber == 1) {
00393 inside = false;
00394 insideskin = false;
00395 tlimit = geombig;
00396 SetParticle( dp->GetDefinition() );
00397 }
00398
00399 theLambdaTable = theTable;
00400 couple = track.GetMaterialCutsCouple();
00401 currentMaterialIndex = couple->GetIndex();
00402 currentKinEnergy = dp->GetKineticEnergy();
00403 currentRange =
00404 theManager->GetRangeFromRestricteDEDX(particle,currentKinEnergy,couple);
00405 lambda0 = GetLambda(currentKinEnergy);
00406
00407
00408 if(inside) return tPathLength;
00409
00410 if(tPathLength > currentRange) tPathLength = currentRange;
00411
00412 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
00413 presafety = sp->GetSafety();
00414
00415
00416
00417
00418
00419
00420 if(currentRange < presafety)
00421 {
00422 inside = true;
00423 return tPathLength;
00424 }
00425
00426 G4StepStatus stepStatus = sp->GetStepStatus();
00427
00428
00429
00430 if (steppingAlgorithm == fUseDistanceToBoundary)
00431 {
00432
00433 GeomLimit(track);
00434
00435
00436 if(currentRange <= presafety)
00437 {
00438 inside = true;
00439 return tPathLength;
00440 }
00441
00442 smallstep += 1.;
00443 insideskin = false;
00444
00445 if((stepStatus == fGeomBoundary) || (stepNumber == 1))
00446 {
00447
00448 if(stepNumber == 1) smallstep = 1.e10;
00449 else smallstep = 1.;
00450
00451
00452
00453 G4double facr = facrange;
00454 if(lambda0 > lambdalimit)
00455 facr *= frscaling1+frscaling2*lambda0/lambdalimit;
00456
00457
00458 if (currentRange > lambda0) tlimit = facr*currentRange;
00459 else tlimit = facr*lambda0;
00460
00461
00462 G4double tgeom = geombig;
00463 if(geomlimit > geommin)
00464 {
00465 if(stepStatus == fGeomBoundary)
00466 tgeom = geomlimit/facgeom;
00467 else
00468 tgeom = 2.*geomlimit/facgeom;
00469 }
00470
00471
00472
00473 G4double rat = currentKinEnergy/MeV ;
00474 rat = 1.e-3/(rat*(10.+rat)) ;
00475
00476 stepmin = rat*lambda0;
00477 skindepth = skin*stepmin;
00478
00479
00480 tlimitmin = lambda0/nstepmax;
00481 if(tlimitmin < stepmin) tlimitmin = 1.01*stepmin;
00482 if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
00483
00484
00485 if(tlimit < tlimitmin) tlimit = tlimitmin;
00486
00487
00488 if(tlimit > tgeom) tlimit = tgeom;
00489 }
00490
00491
00492 if(tlimit < facsafety*presafety) tlimit = facsafety*presafety ;
00493
00494
00495
00496
00497
00498 if((tPathLength < tlimit) && (tPathLength < presafety))
00499 return tPathLength;
00500
00501 G4double tnow = tlimit;
00502
00503 if(geomlimit < geombig) tnow = max(tlimit,facsafety*geomlimit);
00504
00505
00506 if(smallstep < skin)
00507 {
00508 tnow = stepmin;
00509 insideskin = true;
00510 }
00511 else if(geomlimit < geombig)
00512 {
00513 if(geomlimit > skindepth)
00514 {
00515 if(tnow > geomlimit-0.999*skindepth)
00516 tnow = geomlimit-0.999*skindepth;
00517 }
00518 else
00519 {
00520 insideskin = true;
00521 if(tnow > stepmin) tnow = stepmin;
00522 }
00523 }
00524
00525 if(tnow < stepmin) tnow = stepmin;
00526
00527 if(tPathLength > tnow) tPathLength = tnow ;
00528 }
00529
00530
00531 else if(steppingAlgorithm == fUseSafety)
00532 {
00533
00534
00535 if((stepStatus != fGeomBoundary) && (presafety < tlimitminfix))
00536 presafety = safetyHelper->ComputeSafety(sp->GetPosition());
00537
00538
00539 if(currentRange < presafety)
00540 {
00541 inside = true;
00542 return tPathLength;
00543 }
00544
00545 if((stepStatus == fGeomBoundary) || (stepNumber == 1))
00546 {
00547
00548
00549 G4double facr = facrange;
00550 if(lambda0 > lambdalimit)
00551 facr *= frscaling1+frscaling2*lambda0/lambdalimit;
00552
00553
00554 if (currentRange > lambda0) tlimit = facr*currentRange;
00555 else tlimit = facr*lambda0;
00556
00557
00558 tlimitmin = std::max(tlimitminfix,lambda0/nstepmax);
00559 if(tlimit < tlimitmin) tlimit = tlimitmin;
00560 }
00561
00562
00563 if(tlimit < facsafety*presafety) tlimit = facsafety*presafety ;
00564
00565 if(tPathLength > tlimit) tPathLength = tlimit;
00566 }
00567
00568
00569 else
00570 {
00571 if (stepStatus == fGeomBoundary)
00572 {
00573 if (currentRange > lambda0) tlimit = facrange*currentRange;
00574 else tlimit = facrange*lambda0;
00575
00576 if(tlimit < tlimitmin) tlimit = tlimitmin;
00577 if(tPathLength > tlimit) tPathLength = tlimit;
00578 }
00579 }
00580
00581
00582
00583 return tPathLength ;
00584 }
00585
00586
00587
00588 void G4UrbanMscModel90::GeomLimit(const G4Track& track)
00589 {
00590 geomlimit = geombig;
00591
00592
00593 if((track.GetVolume() != 0) &&
00594 (track.GetVolume() != safetyHelper->GetWorldVolume()))
00595 {
00596 G4double cstep = currentRange;
00597 geomlimit = safetyHelper->CheckNextStep(
00598 track.GetStep()->GetPreStepPoint()->GetPosition(),
00599 track.GetMomentumDirection(),
00600 cstep,
00601 presafety);
00602
00603
00604 }
00605 }
00606
00607
00608
00609 G4double G4UrbanMscModel90::ComputeGeomPathLength(G4double)
00610 {
00611 lambdaeff = lambda0;
00612 par1 = -1. ;
00613 par2 = par3 = 0. ;
00614
00615
00616 zPathLength = tPathLength;
00617
00618
00619 if(tPathLength < tlimitminfix) return zPathLength;
00620
00621
00622
00623 if(tPathLength > currentRange)
00624 tPathLength = currentRange ;
00625
00626 G4double tau = tPathLength/lambda0 ;
00627
00628 if ((tau <= tausmall) || insideskin) {
00629 zPathLength = tPathLength;
00630 if(zPathLength > lambda0) zPathLength = lambda0;
00631 return zPathLength;
00632 }
00633
00634 G4double zmean = tPathLength;
00635 if (tPathLength < currentRange*dtrl) {
00636 if(tau < taulim) zmean = tPathLength*(1.-0.5*tau) ;
00637 else zmean = lambda0*(1.-exp(-tau));
00638 } else if(currentKinEnergy < mass) {
00639 par1 = 1./currentRange ;
00640 par2 = 1./(par1*lambda0) ;
00641 par3 = 1.+par2 ;
00642 if(tPathLength < currentRange)
00643 zmean = (1.-exp(par3*log(1.-tPathLength/currentRange)))/(par1*par3) ;
00644 else
00645 zmean = 1./(par1*par3) ;
00646 } else {
00647 G4double T1 = theManager->GetEnergy(particle,currentRange-tPathLength,couple);
00648 G4double lambda1 = GetLambda(T1);
00649
00650 par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
00651 par2 = 1./(par1*lambda0) ;
00652 par3 = 1.+par2 ;
00653 zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ;
00654 }
00655
00656 zPathLength = zmean ;
00657
00658
00659 if(samplez)
00660 {
00661 const G4double ztmax = 0.99, onethird = 1./3. ;
00662 G4double zt = zmean/tPathLength ;
00663
00664 if (tPathLength > stepmin && zt < ztmax)
00665 {
00666 G4double u,cz1;
00667 if(zt >= onethird)
00668 {
00669 G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
00670 cz1 = 1.+cz ;
00671 G4double u0 = cz/cz1 ;
00672 G4double grej ;
00673 do {
00674 u = exp(log(G4UniformRand())/cz1) ;
00675 grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
00676 } while (grej < G4UniformRand()) ;
00677 }
00678 else
00679 {
00680 cz1 = 1./zt-1.;
00681 u = 1.-exp(log(G4UniformRand())/cz1) ;
00682 }
00683 zPathLength = tPathLength*u ;
00684 }
00685 }
00686
00687 if(zPathLength > lambda0) zPathLength = lambda0;
00688
00689 return zPathLength;
00690 }
00691
00692
00693
00694 G4double G4UrbanMscModel90::ComputeTrueStepLength(G4double geomStepLength)
00695 {
00696
00697 if(geomStepLength == zPathLength && tPathLength <= currentRange)
00698 return tPathLength;
00699
00700
00701 zPathLength = geomStepLength;
00702 tPathLength = geomStepLength;
00703 if(geomStepLength < tlimitminfix) return tPathLength;
00704
00705
00706 if((geomStepLength > lambda0*tausmall) && !insideskin)
00707 {
00708 if(par1 < 0.)
00709 tPathLength = -lambda0*log(1.-geomStepLength/lambda0) ;
00710 else
00711 {
00712 if(par1*par3*geomStepLength < 1.)
00713 tPathLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
00714 else
00715 tPathLength = currentRange;
00716 }
00717 }
00718 if(tPathLength < geomStepLength) tPathLength = geomStepLength;
00719
00720 return tPathLength;
00721 }
00722
00723
00724
00725 G4double G4UrbanMscModel90::ComputeTheta0(G4double trueStepLength,
00726 G4double KineticEnergy)
00727 {
00728
00729
00730
00731 const G4double c_highland = 13.6*MeV ;
00732 G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
00733 KineticEnergy*(KineticEnergy+2.*mass)/
00734 ((currentKinEnergy+mass)*(KineticEnergy+mass)));
00735 G4double y = trueStepLength/currentRadLength;
00736 G4double theta0 = c_highland*std::abs(charge)*sqrt(y)/betacp;
00737 y = log(y);
00738 theta0 *= sqrt(1.+y*(0.105+0.0035*y));
00739
00740
00741
00742
00743 theta0 *= 1.-0.24/(Zeff*(Zeff+1.));
00744
00745 return theta0;
00746
00747 }
00748
00749
00750
00751 void G4UrbanMscModel90::SampleScattering(const G4DynamicParticle* dynParticle,
00752 G4double safety)
00753
00754
00755
00756 {
00757 G4double kineticEnergy = dynParticle->GetKineticEnergy();
00758 if((kineticEnergy <= 0.0) || (tPathLength <= tlimitminfix)) return;
00759
00760 G4double cth = SampleCosineTheta(tPathLength,kineticEnergy);
00761
00762 if(cth > 1.) cth = 1.;
00763 if(cth < -1.) cth = -1.;
00764
00765 G4double sth = sqrt((1.0 - cth)*(1.0 + cth));
00766 G4double phi = twopi*G4UniformRand();
00767 G4double dirx = sth*cos(phi);
00768 G4double diry = sth*sin(phi);
00769
00770 G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
00771 G4ThreeVector newDirection(dirx,diry,cth);
00772 newDirection.rotateUz(oldDirection);
00773 fParticleChange->ProposeMomentumDirection(newDirection);
00774
00775 if (latDisplasment && safety > tlimitminfix) {
00776
00777 G4double r = SampleDisplacement();
00778
00779
00780
00781
00782
00783
00784
00785 if(r > 0.)
00786 {
00787 G4double latcorr = LatCorrelation();
00788 if(latcorr > r) latcorr = r;
00789
00790
00791
00792 G4double Phi = 0.;
00793 if(std::abs(r*sth) < latcorr)
00794 Phi = twopi*G4UniformRand();
00795 else
00796 Phi = phi-std::acos(latcorr/(r*sth));
00797 if(Phi < 0.) Phi += twopi;
00798
00799 dirx = std::cos(Phi);
00800 diry = std::sin(Phi);
00801
00802 G4ThreeVector latDirection(dirx,diry,0.0);
00803 latDirection.rotateUz(oldDirection);
00804
00805 G4ThreeVector Position = *(fParticleChange->GetProposedPosition());
00806 G4double fac = 1.;
00807 if(r > safety) {
00808
00809 G4double newsafety = safetyHelper->ComputeSafety(Position);
00810
00811 if(r > newsafety)
00812 fac = newsafety/r ;
00813 }
00814
00815 if(fac > 0.)
00816 {
00817
00818 G4ThreeVector newPosition = Position+fac*r*latDirection;
00819
00820
00821 if(1. == fac) {
00822
00823 safetyHelper->ReLocateWithinVolume(newPosition);
00824
00825
00826 } else {
00827
00828 G4double postsafety = safetyHelper->ComputeSafety(newPosition);
00829
00830
00831 if(postsafety <= 0.) {
00832 safetyHelper->Locate(newPosition, newDirection);
00833
00834
00835 } else {
00836 safetyHelper->ReLocateWithinVolume(newPosition);
00837 }
00838 }
00839 fParticleChange->ProposePosition(newPosition);
00840 }
00841 }
00842 }
00843 }
00844
00845
00846
00847 G4double G4UrbanMscModel90::SampleCosineTheta(G4double trueStepLength,
00848 G4double KineticEnergy)
00849 {
00850 G4double cth = 1. ;
00851 G4double tau = trueStepLength/lambda0 ;
00852
00853 Zeff = couple->GetMaterial()->GetTotNbOfElectPerVolume()/
00854 couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
00855
00856 if(insideskin)
00857 {
00858
00859 G4double mean = trueStepLength/stepmin ;
00860
00861 G4int n = G4Poisson(mean);
00862 if(n > 0)
00863 {
00864 G4double tm = KineticEnergy/electron_mass_c2;
00865
00866 G4double ascr = exp(log(Zeff)/3.)/(137.*sqrt(tm*(tm+2.)));
00867 G4double ascr1 = 1.+0.5*ascr*ascr;
00868 G4double bp1=ascr1+1.;
00869 G4double bm1=ascr1-1.;
00870
00871 G4double ct,st,phi;
00872 G4double sx=0.,sy=0.,sz=0.;
00873 for(G4int i=1; i<=n; i++)
00874 {
00875 ct = ascr1-bp1*bm1/(2.*G4UniformRand()+bm1);
00876 if(ct < -1.) ct = -1.;
00877 if(ct > 1.) ct = 1.;
00878 st = sqrt(1.-ct*ct);
00879 phi = twopi*G4UniformRand();
00880 sx += st*cos(phi);
00881 sy += st*sin(phi);
00882 sz += ct;
00883 }
00884 cth = sz/sqrt(sx*sx+sy*sy+sz*sz);
00885 }
00886 }
00887 else
00888 {
00889 if(trueStepLength >= currentRange*dtrl)
00890 if(par1*trueStepLength < 1.)
00891 tau = -par2*log(1.-par1*trueStepLength) ;
00892
00893
00894 else if(1.-KineticEnergy/currentKinEnergy > taulim)
00895 tau = taubig ;
00896
00897 currentTau = tau ;
00898 lambdaeff = trueStepLength/currentTau;
00899 currentRadLength = couple->GetMaterial()->GetRadlen();
00900
00901 if (tau >= taubig) cth = -1.+2.*G4UniformRand();
00902 else if (tau >= tausmall)
00903 {
00904 G4double b,bx,b1,ebx,eb1;
00905 G4double prob = 0., qprob = 1. ;
00906 G4double a = 1., ea = 0., eaa = 1.;
00907 G4double xmean1 = 1., xmean2 = 0.;
00908 G4double xsi = 3.;
00909
00910 G4double theta0 = ComputeTheta0(trueStepLength,KineticEnergy);
00911
00912
00913 if(theta0 < tausmall) return cth;
00914
00915 G4double sth = sin(0.5*theta0);
00916 a = 0.25/(sth*sth);
00917
00918 G4double xmeanth = exp(-tau);
00919
00920 G4double c = 3. ;
00921 G4double c1 = c-1.;
00922
00923 G4double x0 = 1.-xsi/a ;
00924 if(x0 < 0.)
00925 {
00926
00927 b = exp(tau);
00928 bx = b-1.;
00929 b1 = b+1.;
00930 ebx=exp((c1)*log(bx)) ;
00931 eb1=exp((c1)*log(b1)) ;
00932 }
00933 else
00934 {
00935
00936
00937 c = 2.40-0.027*exp(2.*log(Zeff)/3.);
00938
00939 if(c == 2.) c = 2.+taulim ;
00940 if(c <= 1.) c = 1.+taulim ;
00941 c1 = c-1.;
00942
00943 ea = exp(-xsi) ;
00944 eaa = 1.-ea ;
00945 xmean1 = 1.-(1.-(1.+xsi)*ea)/(eaa*a) ;
00946
00947
00948 b = 1.+(c-xsi)/a ;
00949
00950 b1 = b+1. ;
00951 bx = c/a ;
00952 eb1=exp((c1)*log(b1)) ;
00953 ebx=exp((c1)*log(bx)) ;
00954 xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/(c-2.))/(eb1-ebx) ;
00955
00956 G4double f1x0 = a*ea/eaa ;
00957 G4double f2x0 = c1*eb1*ebx/(eb1-ebx)/exp(c*log(bx)) ;
00958
00959
00960 prob = f2x0/(f1x0+f2x0) ;
00961
00962
00963 qprob = (f1x0+f2x0)*xmeanth/(f2x0*xmean1+f1x0*xmean2) ;
00964 }
00965
00966
00967 if (G4UniformRand() < qprob)
00968 {
00969 if (G4UniformRand() < prob)
00970 cth = 1.+log(ea+G4UniformRand()*eaa)/a ;
00971 else
00972 cth = b-b1*bx/exp(log(ebx-G4UniformRand()*(ebx-eb1))/c1) ;
00973 }
00974 else
00975 {
00976 cth = -1.+2.*G4UniformRand();
00977 }
00978 }
00979 }
00980
00981 return cth ;
00982 }
00983
00984
00985 G4double G4UrbanMscModel90::SampleDisplacement()
00986 {
00987 const G4double kappa = 2.5;
00988 const G4double kappapl1 = kappa+1.;
00989 const G4double kappami1 = kappa-1.;
00990 G4double rmean = 0.0;
00991 if ((currentTau >= tausmall) && !insideskin) {
00992 if (currentTau < taulim) {
00993 rmean = kappa*currentTau*currentTau*currentTau*
00994 (1.-kappapl1*currentTau*0.25)/6. ;
00995
00996 } else {
00997 G4double etau = 0.0;
00998 if (currentTau<taubig) etau = exp(-currentTau);
00999 rmean = -kappa*currentTau;
01000 rmean = -exp(rmean)/(kappa*kappami1);
01001 rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
01002 }
01003 if (rmean>0.) rmean = 2.*lambdaeff*sqrt(rmean/3.0);
01004 else rmean = 0.;
01005 }
01006
01007
01008 if(rmean > 0.) {
01009 G4double zt = (tPathLength-zPathLength)*(tPathLength+zPathLength);
01010 if(zt <= 0.)
01011 rmean = 0.;
01012 else if(rmean*rmean > zt)
01013 rmean = sqrt(zt);
01014 }
01015 return rmean;
01016 }
01017
01018
01019
01020 G4double G4UrbanMscModel90::LatCorrelation()
01021 {
01022 const G4double kappa = 2.5;
01023 const G4double kappami1 = kappa-1.;
01024
01025 G4double latcorr = 0.;
01026 if((currentTau >= tausmall) && !insideskin)
01027 {
01028 if(currentTau < taulim)
01029 latcorr = lambdaeff*kappa*currentTau*currentTau*
01030 (1.-(kappa+1.)*currentTau/3.)/3.;
01031 else
01032 {
01033 G4double etau = 0.;
01034 if(currentTau < taubig) etau = exp(-currentTau);
01035 latcorr = -kappa*currentTau;
01036 latcorr = exp(latcorr)/kappami1;
01037 latcorr += 1.-kappa*etau/kappami1 ;
01038 latcorr *= 2.*lambdaeff/3. ;
01039 }
01040 }
01041
01042 return latcorr;
01043 }
01044
01045
01046
01047 void G4UrbanMscModel90::SampleSecondaries(std::vector<G4DynamicParticle*>*,
01048 const G4MaterialCutsCouple*,
01049 const G4DynamicParticle*,
01050 G4double,
01051 G4double)
01052 {}
01053
01054