G4UrbanMscModel90.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4UrbanMscModel90.cc 66592 2012-12-23 09:34:55Z vnivanch $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:   G4UrbanMscModel90
00034 //
00035 // Author:        V.Ivanchenko clone Laszlo Urban model
00036 //
00037 // Creation date: 07.12.2007
00038 //
00039 // Modifications:
00040 //
00041 //
00042 
00043 // Class Description:
00044 //
00045 // Implementation of the model of multiple scattering based on
00046 // H.W.Lewis Phys Rev 78 (1950) 526 and others
00047 
00048 // -------------------------------------------------------------------
00049 //
00050 
00051 
00052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00054 
00055 #include "G4UrbanMscModel90.hh"
00056 #include "G4PhysicalConstants.hh"
00057 #include "G4SystemOfUnits.hh"
00058 #include "Randomize.hh"
00059 #include "G4Electron.hh"
00060 
00061 #include "G4LossTableManager.hh"
00062 #include "G4ParticleChangeForMSC.hh"
00063 
00064 #include "G4Poisson.hh"
00065 
00066 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00067 
00068 using namespace std;
00069 
00070 G4UrbanMscModel90::G4UrbanMscModel90(const G4String& nam)
00071   : G4VMscModel(nam)
00072 {
00073   taubig        = 8.0;
00074   tausmall      = 1.e-20;
00075   taulim        = 1.e-6;
00076   currentTau    = taulim;
00077   tlimitminfix  = 1.e-6*mm;            
00078   stepmin       = tlimitminfix;
00079   smallstep     = 1.e10;
00080   currentRange  = 0. ;
00081   frscaling2    = 0.25;
00082   frscaling1    = 1.-frscaling2;
00083   tlimit        = 1.e10*mm;
00084   tlimitmin     = 10.*tlimitminfix;            
00085   nstepmax      = 25.;
00086   geombig       = 1.e50*mm;
00087   geommin       = 1.e-3*mm;
00088   geomlimit     = geombig;
00089   presafety     = 0.*mm;
00090   Zeff          = 1.;
00091   particle      = 0;
00092   theManager    = G4LossTableManager::Instance(); 
00093   firstStep     = true; 
00094   inside        = false;  
00095   insideskin    = false;
00096 
00097   skindepth = skin*stepmin;
00098 
00099   mass = proton_mass_c2;
00100   charge = 1.0;
00101   currentKinEnergy = currentRange = currentRadLength = masslimite = masslimitmu 
00102     = lambda0 = lambdaeff = tPathLength = zPathLength = par1 = par2 = par3 = 0;
00103 
00104   currentMaterialIndex = -1;
00105   fParticleChange = 0;
00106   couple = 0;
00107   G4cout << "### G4UrbanMscModel90 is obsolete and will be removed for "
00108          << "the next Geant4 version" << G4endl;
00109 }
00110 
00111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00112 
00113 G4UrbanMscModel90::~G4UrbanMscModel90()
00114 {}
00115 
00116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00117 
00118 void G4UrbanMscModel90::Initialise(const G4ParticleDefinition* p,
00119                                    const G4DataVector&)
00120 {
00121   skindepth = skin*stepmin;
00122 
00123   // set values of some data members
00124   SetParticle(p);
00125 
00126   if(p->GetPDGMass() < MeV) {
00127     G4cout << "### WARNING: G4UrbanMscModel90 model is used for " 
00128            << p->GetParticleName() << " !!! " << G4endl;
00129     G4cout << "###          This model should be used only for heavy particles" 
00130            << G4endl;
00131   }
00132 
00133 
00134   fParticleChange = GetParticleChangeForMSC(p);
00135 }
00136 
00137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00138 
00139 G4double G4UrbanMscModel90::ComputeCrossSectionPerAtom( 
00140                              const G4ParticleDefinition* part,
00141                                    G4double KineticEnergy,
00142                                    G4double AtomicNumber,G4double,
00143                                    G4double, G4double)
00144 {
00145   const G4double sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
00146   const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
00147                             Bohr_radius*Bohr_radius/(hbarc*hbarc);
00148   const G4double epsmin = 1.e-4 , epsmax = 1.e10;
00149 
00150   const G4double Zdat[15] = { 4.,  6., 13., 20., 26., 29., 32., 38., 47.,
00151                              50., 56., 64., 74., 79., 82. };
00152 
00153   const G4double Tdat[22] = { 100*eV,  200*eV,  400*eV,  700*eV,
00154                                1*keV,   2*keV,   4*keV,   7*keV,
00155                               10*keV,  20*keV,  40*keV,  70*keV,
00156                              100*keV, 200*keV, 400*keV, 700*keV,
00157                                1*MeV,   2*MeV,   4*MeV,   7*MeV,
00158                               10*MeV,  20*MeV};
00159 
00160   // corr. factors for e-/e+ lambda for T <= Tlim
00161           G4double celectron[15][22] =
00162           {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
00163             1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
00164             1.112,1.108,1.100,1.093,1.089,1.087            },
00165            {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
00166             1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
00167             1.109,1.105,1.097,1.090,1.086,1.082            },
00168            {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
00169             1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
00170             1.131,1.124,1.113,1.104,1.099,1.098            },
00171            {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
00172             1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
00173             1.112,1.105,1.096,1.089,1.085,1.098            },
00174            {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
00175             1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
00176             1.073,1.070,1.064,1.059,1.056,1.056            },
00177            {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
00178             1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
00179             1.074,1.070,1.063,1.059,1.056,1.052            },
00180            {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
00181             1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
00182             1.068,1.064,1.059,1.054,1.051,1.050            },
00183            {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
00184             1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
00185             1.039,1.037,1.034,1.031,1.030,1.036            },
00186            {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
00187             1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
00188             1.031,1.028,1.024,1.022,1.021,1.024            },
00189            {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
00190             1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
00191             1.020,1.017,1.015,1.013,1.013,1.020            },
00192            {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
00193             1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
00194             0.995,0.993,0.993,0.993,0.993,1.011            },
00195            {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
00196             1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
00197             0.974,0.972,0.973,0.974,0.975,0.987            },
00198            {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
00199             1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
00200             0.950,0.947,0.949,0.952,0.954,0.963            },
00201            {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
00202             1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
00203             0.941,0.938,0.940,0.944,0.946,0.954            },
00204            {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
00205             1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
00206             0.933,0.930,0.933,0.936,0.939,0.949            }};
00207             
00208            G4double cpositron[15][22] = {
00209            {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
00210             1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
00211             1.131,1.126,1.117,1.108,1.103,1.100            },
00212            {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
00213             1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
00214             1.138,1.132,1.122,1.113,1.108,1.102            },
00215            {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
00216             1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
00217             1.203,1.190,1.173,1.159,1.151,1.145            },
00218            {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
00219             1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
00220             1.225,1.210,1.191,1.175,1.166,1.174            },
00221            {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
00222             1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
00223             1.217,1.203,1.184,1.169,1.160,1.151            },
00224            {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
00225             1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
00226             1.237,1.222,1.201,1.184,1.174,1.159            },
00227            {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
00228             1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
00229             1.252,1.234,1.212,1.194,1.183,1.170            },
00230            {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
00231             2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
00232             1.254,1.237,1.214,1.195,1.185,1.179            },
00233            {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
00234             2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
00235             1.312,1.288,1.258,1.235,1.221,1.205            },
00236            {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
00237             2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
00238             1.320,1.294,1.264,1.240,1.226,1.214            },
00239            {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
00240             2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
00241             1.328,1.302,1.270,1.245,1.231,1.233            },
00242            {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
00243             2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
00244             1.361,1.330,1.294,1.267,1.251,1.239            },
00245            {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
00246             3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
00247             1.409,1.372,1.330,1.298,1.280,1.258            },
00248            {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
00249             3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
00250             1.442,1.400,1.354,1.319,1.299,1.272            },
00251            {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
00252             3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
00253             1.456,1.412,1.364,1.328,1.307,1.282            }};
00254 
00255   //data/corrections for T > Tlim  
00256   G4double Tlim = 10.*MeV;
00257   G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
00258                       ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
00259   G4double bg2lim   = Tlim*(Tlim+2.*electron_mass_c2)/
00260                       (electron_mass_c2*electron_mass_c2);
00261 
00262   G4double sig0[15] = {0.2672*barn,  0.5922*barn, 2.653*barn,  6.235*barn,
00263                       11.69*barn  , 13.24*barn  , 16.12*barn, 23.00*barn ,
00264                       35.13*barn  , 39.95*barn  , 50.85*barn, 67.19*barn ,
00265                       91.15*barn  , 104.4*barn  , 113.1*barn};
00266                                        
00267   G4double hecorr[15] = {120.70, 117.50, 105.00, 92.92, 79.23,  74.510,  68.29,
00268                           57.39,  41.97,  36.14, 24.53, 10.21,  -7.855, -16.84,
00269                          -22.30};
00270 
00271   G4double sigma;
00272   SetParticle(part);
00273 
00274   G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
00275 
00276   // correction if particle .ne. e-/e+
00277   // compute equivalent kinetic energy
00278   // lambda depends on p*beta ....
00279 
00280   G4double eKineticEnergy = KineticEnergy;
00281 
00282   if(mass > electron_mass_c2)
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   // interpolate in AtomicNumber and beta2 
00308   G4double c1,c2,cc1,cc2,corr;
00309 
00310   // get bin number in Z
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     // get bin number in T (beta2)
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     //  calculate betasquare values
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00382 
00383 void G4UrbanMscModel90::StartTracking(G4Track* track)
00384 {
00385   SetParticle(track->GetDynamicParticle()->GetDefinition());
00386   firstStep = true; 
00387   inside = false;
00388   insideskin = false;
00389   tlimit = geombig;
00390 }
00391 
00392 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00393 
00394 G4double G4UrbanMscModel90::ComputeTruePathLengthLimit(
00395                              const G4Track& track,
00396                              G4double& currentMinimalStep)
00397 {
00398   tPathLength = currentMinimalStep;
00399   const G4DynamicParticle* dp = track.GetDynamicParticle();
00400   G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
00401   G4StepStatus stepStatus = sp->GetStepStatus();
00402   couple = track.GetMaterialCutsCouple();
00403   SetCurrentCouple(couple); 
00404   currentMaterialIndex = couple->GetIndex();
00405   currentKinEnergy = dp->GetKineticEnergy();
00406   currentRange =  GetRange(particle,currentKinEnergy,couple);
00407   lambda0 = GetTransportMeanFreePath(particle,currentKinEnergy);
00408 
00409   // stop here if small range particle
00410   if(inside || tPathLength < tlimitminfix) { 
00411     return ConvertTrueToGeom(tPathLength, currentMinimalStep); 
00412   }
00413   
00414   if(tPathLength > currentRange) { tPathLength = currentRange; }
00415 
00416   presafety = sp->GetSafety();
00417   /*
00418   G4cout << "G4UrbanMscModel90::ComputeTruePathLengthLimit tPathLength= " 
00419          <<tPathLength<<" safety= " << presafety
00420        << " range= " <<currentRange<<G4endl;
00421   */
00422   // far from geometry boundary
00423   if(currentRange < presafety)
00424     {
00425       inside = true;
00426       return ConvertTrueToGeom(tPathLength, currentMinimalStep);  
00427     }
00428 
00429   // standard  version
00430   //
00431   if (steppingAlgorithm == fUseDistanceToBoundary)
00432     {
00433       //compute geomlimit and presafety 
00434       geomlimit = ComputeGeomLimit(track, presafety, currentRange);
00435    
00436       // is far from boundary
00437       if(currentRange <= presafety)
00438         {
00439           inside = true;
00440           return ConvertTrueToGeom(tPathLength, currentMinimalStep);   
00441         }
00442 
00443       smallstep += 1.;
00444       insideskin = false;
00445 
00446       if(firstStep || stepStatus == fGeomBoundary)
00447         {
00448 
00449           if(firstStep) smallstep = 1.e10;
00450           else  smallstep = 1.;
00451 
00452           // facrange scaling in lambda 
00453           // not so strong step restriction above lambdalimit
00454           G4double facr = facrange;
00455           if(lambda0 > lambdalimit)
00456             facr *= frscaling1+frscaling2*lambda0/lambdalimit;
00457 
00458           // constraint from the physics
00459           if (currentRange > lambda0) tlimit = facr*currentRange;
00460           else                        tlimit = facr*lambda0;
00461 
00462           // constraint from the geometry (if tlimit above is too big)
00463           G4double tgeom = geombig; 
00464           if(geomlimit > geommin)
00465             {
00466               if(stepStatus == fGeomBoundary)  
00467                 tgeom = geomlimit/facgeom;
00468               else
00469                 tgeom = 2.*geomlimit/facgeom;
00470             }
00471 
00472           //define stepmin here (it depends on lambda!)
00473           //rough estimation of lambda_elastic/lambda_transport
00474           G4double rat = currentKinEnergy/MeV ;
00475           rat = 1.e-3/(rat*(10.+rat)) ;
00476           //stepmin ~ lambda_elastic
00477           stepmin = rat*lambda0;
00478           skindepth = skin*stepmin;
00479 
00480           //define tlimitmin
00481           tlimitmin = lambda0/nstepmax;
00482           if(tlimitmin < stepmin) tlimitmin = 1.01*stepmin;
00483           if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
00484 
00485           //lower limit for tlimit
00486           if(tlimit < tlimitmin) tlimit = tlimitmin;
00487 
00488           //check against geometry limit
00489           if(tlimit > tgeom) tlimit = tgeom;
00490         }
00491 
00492       //if track starts far from boundaries increase tlimit!
00493       if(tlimit < facsafety*presafety) tlimit = facsafety*presafety ;
00494 
00495       //  G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit  
00496       //     << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
00497 
00498       // shortcut
00499       if((tPathLength < tlimit) && (tPathLength < presafety))
00500         return ConvertTrueToGeom(tPathLength, currentMinimalStep);   
00501 
00502       G4double tnow = tlimit;
00503       // optimization ...
00504       if(geomlimit < geombig) tnow = max(tlimit,facsafety*geomlimit);
00505    
00506       // step reduction near to boundary
00507       if(smallstep < skin)
00508         {
00509           tnow = stepmin;
00510           insideskin = true;
00511         }
00512       else if(geomlimit < geombig)
00513         {
00514           if(geomlimit > skindepth)
00515             {
00516               if(tnow > geomlimit-0.999*skindepth)
00517                 tnow = geomlimit-0.999*skindepth;
00518             }
00519           else
00520             {
00521               insideskin = true;
00522               if(tnow > stepmin) tnow = stepmin;
00523             }
00524         }
00525 
00526       if(tnow < stepmin) tnow = stepmin;
00527 
00528       if(tPathLength > tnow) tPathLength = tnow ; 
00529     }
00530     // for 'normal' simulation with or without magnetic field 
00531     //  there no small step/single scattering at boundaries
00532   else if(steppingAlgorithm == fUseSafety)
00533     {
00534       // compute presafety again if presafety <= 0 and no boundary
00535       // i.e. when it is needed for optimization purposes
00536       if((stepStatus != fGeomBoundary) && (presafety < tlimitminfix)) 
00537         presafety = ComputeSafety(sp->GetPosition(),tPathLength); 
00538 
00539       // is far from boundary
00540       if(currentRange < presafety)
00541         {
00542           inside = true;
00543           return ConvertTrueToGeom(tPathLength, currentMinimalStep);  
00544         }
00545 
00546       if(firstStep || stepStatus == fGeomBoundary)
00547         { 
00548           // facrange scaling in lambda 
00549           // not so strong step restriction above lambdalimit
00550           G4double facr = facrange;
00551           if(lambda0 > lambdalimit)
00552             facr *= frscaling1+frscaling2*lambda0/lambdalimit;
00553 
00554           // constraint from the physics
00555           if (currentRange > lambda0) tlimit = facr*currentRange;
00556           else                        tlimit = facr*lambda0;
00557 
00558           //lower limit for tlimit
00559           tlimitmin = std::max(tlimitminfix,lambda0/nstepmax);
00560           if(tlimit < tlimitmin) tlimit = tlimitmin;
00561         }
00562 
00563       //if track starts far from boundaries increase tlimit!
00564       if(tlimit < facsafety*presafety) tlimit = facsafety*presafety ;
00565 
00566       if(tPathLength > tlimit) tPathLength = tlimit;
00567     }
00568   
00569   // version similar to 7.1 (needed for some experiments)
00570   else
00571     {
00572       if (stepStatus == fGeomBoundary)
00573         {
00574           if (currentRange > lambda0) tlimit = facrange*currentRange;
00575           else                        tlimit = facrange*lambda0;
00576 
00577           if(tlimit < tlimitmin) tlimit = tlimitmin;
00578           if(tPathLength > tlimit) tPathLength = tlimit;
00579         }
00580     }
00581   //  G4cout << "tPathLength= " << tPathLength << "  geomlimit= " << geomlimit 
00582   //     << " currentMinimalStep= " << currentMinimalStep << G4endl;
00583   return ConvertTrueToGeom(tPathLength, currentMinimalStep);
00584 }
00585 
00586 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00587 
00588 G4double G4UrbanMscModel90::ComputeGeomPathLength(G4double)
00589 {
00590   firstStep = false; 
00591   lambdaeff = lambda0;
00592   par1 = -1. ;  
00593   par2 = par3 = 0. ;  
00594 
00595   //  do the true -> geom transformation
00596   zPathLength = tPathLength;
00597 
00598   // z = t for very small tPathLength
00599   if(tPathLength < tlimitminfix) return zPathLength;
00600 
00601   // this correction needed to run MSC with eIoni and eBrem inactivated
00602   // and makes no harm for a normal run
00603   if(tPathLength > currentRange)
00604     tPathLength = currentRange ;
00605 
00606   G4double tau   = tPathLength/lambda0 ;
00607 
00608   if ((tau <= tausmall) || insideskin) {
00609     zPathLength  = tPathLength;
00610     if(zPathLength > lambda0) zPathLength = lambda0;
00611     return zPathLength;
00612   }
00613 
00614   G4double zmean = tPathLength;
00615   if (tPathLength < currentRange*dtrl) {
00616     if(tau < taulim) zmean = tPathLength*(1.-0.5*tau) ;
00617     else             zmean = lambda0*(1.-exp(-tau));
00618   } else if(currentKinEnergy < mass || tPathLength == currentRange) {
00619     par1 = 1./currentRange ;
00620     par2 = 1./(par1*lambda0) ;
00621     par3 = 1.+par2 ;
00622     if(tPathLength < currentRange)
00623       zmean = (1.-exp(par3*log(1.-tPathLength/currentRange)))/(par1*par3) ;
00624     else
00625       zmean = 1./(par1*par3) ;
00626   } else {
00627     G4double T1 = GetEnergy(particle,currentRange-tPathLength,couple);
00628     G4double lambda1 = GetTransportMeanFreePath(particle,T1);
00629 
00630     par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
00631     par2 = 1./(par1*lambda0) ;
00632     par3 = 1.+par2 ;
00633     zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ;
00634   }
00635 
00636   zPathLength = zmean ;
00637 
00638   //  sample z
00639   if(samplez)
00640   {
00641     const G4double  ztmax = 0.99, onethird = 1./3. ;
00642     G4double zt = zmean/tPathLength ;
00643 
00644     if (tPathLength > stepmin && zt < ztmax)              
00645     {
00646       G4double u,cz1;
00647       if(zt >= onethird)
00648       {
00649         G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
00650         cz1 = 1.+cz ;
00651         G4double u0 = cz/cz1 ;
00652         G4double grej ;
00653         do {
00654             u = exp(log(G4UniformRand())/cz1) ;
00655             grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
00656            } while (grej < G4UniformRand()) ;
00657       }
00658       else
00659       {
00660         cz1 = 1./zt-1.;
00661         u = 1.-exp(log(G4UniformRand())/cz1) ;
00662       }
00663       zPathLength = tPathLength*u ;
00664     }
00665   }
00666 
00667   if(zPathLength > lambda0) zPathLength = lambda0;
00668 
00669   return zPathLength;
00670 }
00671 
00672 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00673 
00674 G4double G4UrbanMscModel90::ComputeTrueStepLength(G4double geomStepLength)
00675 {
00676   // step defined other than transportation 
00677   if(geomStepLength == zPathLength && tPathLength <= currentRange)
00678     { return tPathLength; }
00679 
00680   // t = z for very small step
00681   zPathLength = geomStepLength;
00682   tPathLength = geomStepLength;
00683   if(geomStepLength < tlimitminfix) return tPathLength;
00684   
00685   // recalculation
00686   if((geomStepLength > lambda0*tausmall) && !insideskin)
00687   {
00688     if(par1 <  0.)
00689       tPathLength = -lambda0*log(1.-geomStepLength/lambda0) ;
00690     else 
00691     {
00692       if(par1*par3*geomStepLength < 1.)
00693         tPathLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
00694       else 
00695         tPathLength = currentRange;
00696     }  
00697   }
00698   if(tPathLength < geomStepLength) tPathLength = geomStepLength;
00699 
00700   return tPathLength;
00701 }
00702 
00703 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00704 
00705 G4double G4UrbanMscModel90::ComputeTheta0(G4double trueStepLength,
00706                                         G4double KineticEnergy)
00707 {
00708   // for all particles take the width of the central part
00709   //  from a  parametrization similar to the Highland formula
00710   // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
00711   const G4double c_highland = 13.6*MeV ;
00712   G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
00713                          KineticEnergy*(KineticEnergy+2.*mass)/
00714                       ((currentKinEnergy+mass)*(KineticEnergy+mass)));
00715   G4double y = trueStepLength/currentRadLength;
00716   G4double theta0 = c_highland*std::abs(charge)*sqrt(y)/betacp;
00717            y = log(y);
00718            theta0 *= sqrt(1.+y*(0.105+0.0035*y));
00719 
00720   //correction for small Zeff (based on high energy
00721   // proton scattering  data)
00722   // see G.Shen at al. Phys.Rev.D20(1979) p.1584
00723   theta0 *= 1.-0.24/(Zeff*(Zeff+1.));
00724 
00725   return theta0;
00726 
00727 }
00728 
00729 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00730 
00731 G4ThreeVector& 
00732 G4UrbanMscModel90::SampleScattering(const G4ThreeVector& oldDirection,
00733                                     G4double safety)
00734 {
00735   fDisplacement.set(0.0,0.0,0.0);
00736   G4double kineticEnergy = currentKinEnergy;
00737   if (tPathLength > currentRange*dtrl) {
00738     kineticEnergy = GetEnergy(particle,currentRange-tPathLength,couple);
00739   } else {
00740     kineticEnergy -= tPathLength*GetDEDX(particle,currentKinEnergy,couple);
00741   }
00742   if((kineticEnergy <= eV) || (tPathLength <= tlimitminfix) ||
00743      (tPathLength/tausmall < lambda0) ) { return fDisplacement; }
00744 
00745   G4double cth  = SampleCosineTheta(tPathLength,kineticEnergy);
00746   // protection against 'bad' cth values
00747   if(std::fabs(cth) > 1.) { return fDisplacement; }
00748 
00749   const G4double checkEnergy = GeV;
00750   if(kineticEnergy > checkEnergy && cth < 0.0 
00751      && tPathLength < taulim*lambda0) {
00752     G4ExceptionDescription ed;
00753     ed << particle->GetParticleName()
00754        << " E(MeV)= " << kineticEnergy/MeV
00755        << " Step(mm)= " << tPathLength/mm
00756        << " in " << CurrentCouple()->GetMaterial()->GetName()
00757        << " CosTheta= " << cth 
00758        << " is too big - the angle is set to zero" << G4endl;
00759     G4Exception("G4UrbanMscModel90::SampleScattering","em0004",JustWarning,
00760                 ed,"");
00761     return fDisplacement;
00762   }
00763 
00764   G4double sth  = sqrt((1.0 - cth)*(1.0 + cth));
00765   G4double phi  = twopi*G4UniformRand();
00766   G4double dirx = sth*cos(phi);
00767   G4double diry = sth*sin(phi);
00768 
00769   // G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
00770   G4ThreeVector newDirection(dirx,diry,cth);
00771   newDirection.rotateUz(oldDirection);
00772   fParticleChange->ProposeMomentumDirection(newDirection);
00773 
00774   if (latDisplasment && safety > tlimitminfix) {
00775 
00776     G4double r = SampleDisplacement();
00777 /*
00778     G4cout << "G4UrbanMscModel90::SampleSecondaries: e(MeV)= " << kineticEnergy
00779            << " sinTheta= " << sth << " r(mm)= " << r
00780            << " trueStep(mm)= " << truestep 
00781            << " geomStep(mm)= " << zPathLength
00782            << G4endl;
00783 */
00784     if(r > 0.)
00785       {
00786         G4double latcorr = LatCorrelation();
00787         if(latcorr > r) latcorr = r;
00788 
00789         // sample direction of lateral displacement
00790         // compute it from the lateral correlation
00791         G4double Phi = 0.;
00792         if(std::abs(r*sth) < latcorr) {
00793           Phi  = twopi*G4UniformRand();
00794         } else {
00795           G4double psi = std::acos(latcorr/(r*sth));
00796           if(G4UniformRand() < 0.5) Phi = phi+psi;
00797           else                      Phi = phi-psi;
00798         }
00799 
00800         dirx = std::cos(Phi);
00801         diry = std::sin(Phi);
00802 
00803         fDisplacement.set(r*dirx,r*diry,0.0);
00804         fDisplacement.rotateUz(oldDirection);
00805       }
00806   }
00807   return fDisplacement;
00808 }
00809 
00810 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00811 
00812 G4double G4UrbanMscModel90::SampleCosineTheta(G4double trueStepLength,
00813                                             G4double KineticEnergy)
00814 {
00815   G4double cth = 1. ;
00816   G4double tau = trueStepLength/lambda0 ;
00817 
00818   Zeff = couple->GetMaterial()->GetTotNbOfElectPerVolume()/
00819          couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
00820 
00821   if(insideskin)
00822   {
00823     //no scattering, single or plural scattering
00824     G4double mean = trueStepLength/stepmin ;
00825 
00826     G4int n = G4Poisson(mean);
00827     if(n > 0)
00828     {
00829       G4double tm = KineticEnergy/electron_mass_c2;
00830       // ascr - screening parameter
00831       G4double ascr = exp(log(Zeff)/3.)/(137.*sqrt(tm*(tm+2.)));
00832       G4double ascr1 = 1.+0.5*ascr*ascr;
00833       G4double bp1=ascr1+1.;
00834       G4double bm1=ascr1-1.;
00835       // single scattering from screened Rutherford x-section
00836       G4double ct,st,phi;
00837       G4double sx=0.,sy=0.,sz=0.;
00838       for(G4int i=1; i<=n; i++)
00839       {
00840         ct = ascr1-bp1*bm1/(2.*G4UniformRand()+bm1);
00841         if(ct < -1.) ct = -1.;
00842         if(ct >  1.) ct =  1.; 
00843         st = sqrt(1.-ct*ct);
00844         phi = twopi*G4UniformRand();
00845         sx += st*cos(phi);
00846         sy += st*sin(phi);
00847         sz += ct;
00848       }
00849         cth = sz/sqrt(sx*sx+sy*sy+sz*sz);
00850     }
00851   }
00852   else
00853   {
00854     if(trueStepLength >= currentRange*dtrl) {
00855       if(par1*trueStepLength < 1.)
00856         tau = -par2*log(1.-par1*trueStepLength) ;
00857       // for the case if ioni/brems are inactivated
00858       // see the corresponding condition in ComputeGeomPathLength 
00859       else if(1.-KineticEnergy/currentKinEnergy > taulim)
00860         tau = taubig ;
00861     }
00862     currentTau = tau ;
00863     lambdaeff = trueStepLength/currentTau;
00864     currentRadLength = couple->GetMaterial()->GetRadlen();
00865 
00866     if (tau >= taubig) cth = -1.+2.*G4UniformRand();
00867     else if (tau >= tausmall)
00868     {
00869       G4double b,bx,b1,ebx,eb1;
00870       G4double prob = 0., qprob = 1. ;
00871       G4double a = 1., ea = 0., eaa = 1.;
00872       G4double xmean1 = 1., xmean2 = 0.;
00873       G4double xsi = 3.;
00874 
00875       G4double theta0 = ComputeTheta0(trueStepLength,KineticEnergy);
00876 
00877       // protexction for very small angles
00878       if(theta0*theta0 < tausmall) return cth;
00879 
00880       G4double sth = sin(0.5*theta0);
00881       a = 0.25/(sth*sth);
00882 
00883       G4double xmeanth = exp(-tau);
00884 
00885       G4double c = 3. ;         
00886       G4double c1 = c-1.;
00887 
00888       G4double x0 = 1.-xsi/a ;
00889       if(x0 < 0.)
00890       {
00891         // 1 model function
00892         b = exp(tau);
00893         bx = b-1.;
00894         b1 = b+1.;
00895         ebx=exp((c1)*log(bx)) ;
00896         eb1=exp((c1)*log(b1)) ;
00897       }
00898       else
00899       {
00900         //empirical tail parameter 
00901         // based some exp. data
00902         c = 2.40-0.027*exp(2.*log(Zeff)/3.);
00903 
00904         if(c == 2.) c = 2.+taulim ;
00905         if(c <= 1.) c = 1.+taulim ;
00906         c1 = c-1.;
00907 
00908         ea = exp(-xsi) ;
00909         eaa = 1.-ea ;
00910         xmean1 = 1.-(1.-(1.+xsi)*ea)/(eaa*a) ; 
00911 
00912         // from the continuity of the 1st derivative at x=x0
00913         b = 1.+(c-xsi)/a ;
00914 
00915         b1 = b+1. ;
00916         bx = c/a ;
00917         eb1=exp((c1)*log(b1)) ;
00918         ebx=exp((c1)*log(bx)) ;
00919         xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/(c-2.))/(eb1-ebx) ;
00920 
00921         G4double f1x0 = a*ea/eaa ;
00922         G4double f2x0 = c1*eb1*ebx/(eb1-ebx)/exp(c*log(bx)) ;
00923 
00924         // from continuity at x=x0
00925         prob = f2x0/(f1x0+f2x0) ;
00926 
00927         // from xmean = xmeanth
00928         qprob = (f1x0+f2x0)*xmeanth/(f2x0*xmean1+f1x0*xmean2) ;
00929       }
00930 
00931       // sampling of costheta
00932       if (G4UniformRand() < qprob)
00933       {
00934         if (G4UniformRand() < prob)
00935            cth = 1.+log(ea+G4UniformRand()*eaa)/a ;
00936         else
00937            cth = b-b1*bx/exp(log(ebx-G4UniformRand()*(ebx-eb1))/c1) ;
00938       }
00939       else
00940       {
00941         cth = -1.+2.*G4UniformRand();
00942       }
00943     }
00944   }  
00945 
00946   return cth ;
00947 }
00948 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00949 
00950 G4double G4UrbanMscModel90::SampleDisplacement()
00951 {
00952   const G4double kappa = 2.5;
00953   const G4double kappapl1 = kappa+1.;
00954   const G4double kappami1 = kappa-1.;
00955   G4double rmean = 0.0;
00956   if ((currentTau >= tausmall) && !insideskin) {
00957     if (currentTau < taulim) {
00958       rmean = kappa*currentTau*currentTau*currentTau*
00959              (1.-kappapl1*currentTau*0.25)/6. ;
00960 
00961     } else {
00962       G4double etau = 0.0;
00963       if (currentTau<taubig) etau = exp(-currentTau);
00964       rmean = -kappa*currentTau;
00965       rmean = -exp(rmean)/(kappa*kappami1);
00966       rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
00967     }
00968     if (rmean>0.) rmean = 2.*lambdaeff*sqrt(rmean/3.0);
00969     else          rmean = 0.;
00970   }
00971 
00972   // protection against z > t ...........................
00973   if(rmean > 0.) {
00974     G4double zt = (tPathLength-zPathLength)*(tPathLength+zPathLength);
00975     if(zt <= 0.)
00976       rmean = 0.;
00977     else if(rmean*rmean > zt)
00978       rmean = sqrt(zt);
00979   }
00980   return rmean;
00981 }
00982 
00983 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00984 
00985 G4double G4UrbanMscModel90::LatCorrelation()
00986 {
00987   const G4double kappa = 2.5;
00988   const G4double kappami1 = kappa-1.;
00989 
00990   G4double latcorr = 0.;
00991   if((currentTau >= tausmall) && !insideskin)
00992   {
00993     if(currentTau < taulim)
00994       latcorr = lambdaeff*kappa*currentTau*currentTau*
00995                 (1.-(kappa+1.)*currentTau/3.)/3.;
00996     else
00997     {
00998       G4double etau = 0.;
00999       if(currentTau < taubig) etau = exp(-currentTau);
01000       latcorr = -kappa*currentTau;
01001       latcorr = exp(latcorr)/kappami1;
01002       latcorr += 1.-kappa*etau/kappami1 ;
01003       latcorr *= 2.*lambdaeff/3. ;
01004     }
01005   }
01006 
01007   return latcorr;
01008 }
01009 
01010 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
01011 
01012 void G4UrbanMscModel90::SampleSecondaries(std::vector<G4DynamicParticle*>*,
01013                                           const G4MaterialCutsCouple*,
01014                                           const G4DynamicParticle*,
01015                                           G4double,
01016                                           G4double)
01017 {}
01018 
01019 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

Generated on Mon May 27 17:50:08 2013 for Geant4 by  doxygen 1.4.7