G4UrbanMscModel96.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:   $
00027 // GEANT4 tag $Name:  $
00028 //
00029 // -------------------------------------------------------------------
00030 //   
00031 // GEANT4 Class file
00032 //    
00033 //
00034 // File name:   G4UrbanMscModel96
00035 //
00036 // Author:      Laszlo Urban
00037 //
00038 // Creation date: 26.09.2012
00039 //
00040 // Created from G4UrbanMscModel95
00041 //
00042 // New parametrization for theta0
00043 // Correction for very small step length
00044 //
00045 // Class Description:
00046 //
00047 // Implementation of the model of multiple scattering based on
00048 // H.W.Lewis Phys Rev 78 (1950) 526 and others
00049 
00050 // -------------------------------------------------------------------
00051 // In its present form the model can be  used for simulation 
00052 //   of the e-/e+ multiple scattering
00053 //
00054 
00055 
00056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00058 
00059 #include "G4UrbanMscModel96.hh"
00060 #include "G4PhysicalConstants.hh"
00061 #include "G4SystemOfUnits.hh"
00062 #include "Randomize.hh"
00063 #include "G4Electron.hh"
00064 #include "G4LossTableManager.hh"
00065 #include "G4ParticleChangeForMSC.hh"
00066 
00067 #include "G4Poisson.hh"
00068 #include "globals.hh"
00069 #include "G4Pow.hh"
00070 
00071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00072 
00073 using namespace std;
00074 
00075 G4UrbanMscModel96::G4UrbanMscModel96(const G4String& nam)
00076   : G4VMscModel(nam)
00077 {
00078   masslimite    = 0.6*MeV;
00079   lambdalimit   = 1.*mm;
00080   fr            = 0.02;
00081   taubig        = 8.0;
00082   tausmall      = 1.e-16;
00083   taulim        = 1.e-6;
00084   currentTau    = taulim;
00085   tlimitminfix  = 1.e-6*mm;            
00086   stepmin       = tlimitminfix;
00087   smallstep     = 1.e10;
00088   currentRange  = 0. ;
00089   rangeinit     = 0.;
00090   tlimit        = 1.e10*mm;
00091   tlimitmin     = 10.*tlimitminfix;            
00092   tgeom         = 1.e50*mm;
00093   geombig       = 1.e50*mm;
00094   geommin       = 1.e-3*mm;
00095   geomlimit     = geombig;
00096   presafety     = 0.*mm;
00097   //facsafety     = 0.50 ;
00098                           
00099   y             = 0.;
00100 
00101   Zold          = 0.;
00102   Zeff          = 1.;
00103   Z2            = 1.;                
00104   Z23           = 1.;                    
00105   lnZ           = 0.;
00106   coeffth1      = 0.;
00107   coeffth2      = 0.;
00108   coeffc1       = 0.;
00109   coeffc2       = 0.;
00110   coeffc3       = 0.;
00111   coeffc4       = 0.;
00112 
00113   theta0max     = pi/6.;
00114   rellossmax    = 0.50;
00115   third         = 1./3.;
00116   particle      = 0;
00117   theManager    = G4LossTableManager::Instance(); 
00118   firstStep     = true; 
00119   inside        = false;  
00120   insideskin    = false;
00121 
00122   skindepth = skin*stepmin;
00123 
00124   mass = proton_mass_c2;
00125   charge = ChargeSquare = 1.0;
00126   currentKinEnergy = currentRadLength = lambda0 = lambdaeff = tPathLength 
00127     = zPathLength = par1 = par2 = par3 = 0;
00128 
00129   currentMaterialIndex = -1;
00130   fParticleChange = 0;
00131   couple = 0;
00132   SetSampleZ(false);
00133 }
00134 
00135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00136 
00137 G4UrbanMscModel96::~G4UrbanMscModel96()
00138 {}
00139 
00140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00141 
00142 void G4UrbanMscModel96::Initialise(const G4ParticleDefinition* p,
00143                                    const G4DataVector&)
00144 {
00145   skindepth = skin*stepmin;
00146   //  trackID = -1;
00147 
00148   // set values of some data members
00149   SetParticle(p);
00150   /*
00151   if(p->GetPDGMass() > MeV) {
00152     G4cout << "### WARNING: G4UrbanMscModel96 model is used for " 
00153            << p->GetParticleName() << " !!! " << G4endl;
00154     G4cout << "###          This model should be used only for e+-" 
00155            << G4endl;
00156   }
00157   */
00158   fParticleChange = GetParticleChangeForMSC(p);
00159 
00160   //samplez = true;
00161   //G4cout << "### G4UrbanMscModel96::Initialise done!" << G4endl;
00162 }
00163 
00164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00165 
00166 G4double G4UrbanMscModel96::ComputeCrossSectionPerAtom( 
00167                              const G4ParticleDefinition* part,
00168                                    G4double KineticEnergy,
00169                                    G4double AtomicNumber,G4double,
00170                                    G4double, G4double)
00171 {
00172   static const G4double sigmafactor = 
00173     twopi*classic_electr_radius*classic_electr_radius;
00174   static const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
00175                             Bohr_radius*Bohr_radius/(hbarc*hbarc);
00176   static const G4double epsmin = 1.e-4 , epsmax = 1.e10;
00177 
00178   static const G4double Zdat[15] = { 4.,  6., 13., 20., 26., 29., 32., 38., 47.,
00179                                      50., 56., 64., 74., 79., 82. };
00180 
00181   static const G4double Tdat[22] = { 100*eV,  200*eV,  400*eV,  700*eV,
00182                                      1*keV,   2*keV,   4*keV,   7*keV,
00183                                      10*keV,  20*keV,  40*keV,  70*keV,
00184                                      100*keV, 200*keV, 400*keV, 700*keV,
00185                                      1*MeV,   2*MeV,   4*MeV,   7*MeV,
00186                                      10*MeV,  20*MeV};
00187 
00188   // corr. factors for e-/e+ lambda for T <= Tlim
00189   static const G4double celectron[15][22] =
00190           {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
00191             1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
00192             1.112,1.108,1.100,1.093,1.089,1.087            },
00193            {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
00194             1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
00195             1.109,1.105,1.097,1.090,1.086,1.082            },
00196            {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
00197             1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
00198             1.131,1.124,1.113,1.104,1.099,1.098            },
00199            {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
00200             1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
00201             1.112,1.105,1.096,1.089,1.085,1.098            },
00202            {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
00203             1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
00204             1.073,1.070,1.064,1.059,1.056,1.056            },
00205            {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
00206             1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
00207             1.074,1.070,1.063,1.059,1.056,1.052            },
00208            {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
00209             1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
00210             1.068,1.064,1.059,1.054,1.051,1.050            },
00211            {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
00212             1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
00213             1.039,1.037,1.034,1.031,1.030,1.036            },
00214            {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
00215             1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
00216             1.031,1.028,1.024,1.022,1.021,1.024            },
00217            {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
00218             1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
00219             1.020,1.017,1.015,1.013,1.013,1.020            },
00220            {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
00221             1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
00222             0.995,0.993,0.993,0.993,0.993,1.011            },
00223            {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
00224             1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
00225             0.974,0.972,0.973,0.974,0.975,0.987            },
00226            {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
00227             1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
00228             0.950,0.947,0.949,0.952,0.954,0.963            },
00229            {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
00230             1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
00231             0.941,0.938,0.940,0.944,0.946,0.954            },
00232            {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
00233             1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
00234             0.933,0.930,0.933,0.936,0.939,0.949            }};
00235             
00236   static const G4double cpositron[15][22] = {
00237            {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
00238             1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
00239             1.131,1.126,1.117,1.108,1.103,1.100            },
00240            {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
00241             1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
00242             1.138,1.132,1.122,1.113,1.108,1.102            },
00243            {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
00244             1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
00245             1.203,1.190,1.173,1.159,1.151,1.145            },
00246            {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
00247             1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
00248             1.225,1.210,1.191,1.175,1.166,1.174            },
00249            {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
00250             1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
00251             1.217,1.203,1.184,1.169,1.160,1.151            },
00252            {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
00253             1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
00254             1.237,1.222,1.201,1.184,1.174,1.159            },
00255            {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
00256             1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
00257             1.252,1.234,1.212,1.194,1.183,1.170            },
00258            {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
00259             2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
00260             1.254,1.237,1.214,1.195,1.185,1.179            },
00261            {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
00262             2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
00263             1.312,1.288,1.258,1.235,1.221,1.205            },
00264            {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
00265             2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
00266             1.320,1.294,1.264,1.240,1.226,1.214            },
00267            {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
00268             2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
00269             1.328,1.302,1.270,1.245,1.231,1.233            },
00270            {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
00271             2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
00272             1.361,1.330,1.294,1.267,1.251,1.239            },
00273            {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
00274             3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
00275             1.409,1.372,1.330,1.298,1.280,1.258            },
00276            {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
00277             3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
00278             1.442,1.400,1.354,1.319,1.299,1.272            },
00279            {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
00280             3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
00281             1.456,1.412,1.364,1.328,1.307,1.282            }};
00282 
00283   //data/corrections for T > Tlim  
00284   static const G4double Tlim = 10.*MeV;
00285   static const G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
00286                       ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
00287   static const G4double bg2lim   = Tlim*(Tlim+2.*electron_mass_c2)/
00288                       (electron_mass_c2*electron_mass_c2);
00289 
00290   static const G4double sig0[15] = {
00291                      0.2672*barn,  0.5922*barn, 2.653*barn,  6.235*barn,
00292                       11.69*barn  , 13.24*barn  , 16.12*barn, 23.00*barn ,
00293                       35.13*barn  , 39.95*barn  , 50.85*barn, 67.19*barn ,
00294                       91.15*barn  , 104.4*barn  , 113.1*barn};
00295                                        
00296   static const G4double hecorr[15] = {
00297                          120.70, 117.50, 105.00, 92.92, 79.23,  74.510,  68.29,
00298                           57.39,  41.97,  36.14, 24.53, 10.21,  -7.855, -16.84,
00299                          -22.30};
00300 
00301   G4double sigma;
00302   SetParticle(part);
00303 
00304   Z23 = G4Pow::GetInstance()->Z23(G4lrint(AtomicNumber));
00305 
00306   // correction if particle .ne. e-/e+
00307   // compute equivalent kinetic energy
00308   // lambda depends on p*beta ....
00309 
00310   G4double eKineticEnergy = KineticEnergy;
00311 
00312   if(mass > electron_mass_c2)
00313   {
00314      G4double TAU = KineticEnergy/mass ;
00315      G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
00316      G4double w = c-2. ;
00317      G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
00318      eKineticEnergy = electron_mass_c2*tau ;
00319   }
00320 
00321   G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
00322   G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00323                                  /(eTotalEnergy*eTotalEnergy);
00324   G4double bg2   = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00325                                  /(electron_mass_c2*electron_mass_c2);
00326 
00327   G4double eps = epsfactor*bg2/Z23;
00328 
00329   if     (eps<epsmin)  sigma = 2.*eps*eps;
00330   else if(eps<epsmax)  sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
00331   else                 sigma = log(2.*eps)-1.+1./eps;
00332 
00333   sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
00334 
00335   // interpolate in AtomicNumber and beta2 
00336   G4double c1,c2,cc1,cc2,corr;
00337 
00338   // get bin number in Z
00339   G4int iZ = 14;
00340   while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
00341   if (iZ==14)                               iZ = 13;
00342   if (iZ==-1)                               iZ = 0 ;
00343 
00344   G4double ZZ1 = Zdat[iZ];
00345   G4double ZZ2 = Zdat[iZ+1];
00346   G4double ratZ = (AtomicNumber-ZZ1)*(AtomicNumber+ZZ1)/
00347                   ((ZZ2-ZZ1)*(ZZ2+ZZ1));
00348 
00349   if(eKineticEnergy <= Tlim) 
00350   {
00351     // get bin number in T (beta2)
00352     G4int iT = 21;
00353     while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
00354     if(iT==21)                                  iT = 20;
00355     if(iT==-1)                                  iT = 0 ;
00356 
00357     //  calculate betasquare values
00358     G4double T = Tdat[iT],   E = T + electron_mass_c2;
00359     G4double b2small = T*(E+electron_mass_c2)/(E*E);
00360 
00361     T = Tdat[iT+1]; E = T + electron_mass_c2;
00362     G4double b2big = T*(E+electron_mass_c2)/(E*E);
00363     G4double ratb2 = (beta2-b2small)/(b2big-b2small);
00364 
00365     if (charge < 0.)
00366     {
00367        c1 = celectron[iZ][iT];
00368        c2 = celectron[iZ+1][iT];
00369        cc1 = c1+ratZ*(c2-c1);
00370 
00371        c1 = celectron[iZ][iT+1];
00372        c2 = celectron[iZ+1][iT+1];
00373        cc2 = c1+ratZ*(c2-c1);
00374 
00375        corr = cc1+ratb2*(cc2-cc1);
00376 
00377        sigma *= sigmafactor/corr;
00378     }
00379     else              
00380     {
00381        c1 = cpositron[iZ][iT];
00382        c2 = cpositron[iZ+1][iT];
00383        cc1 = c1+ratZ*(c2-c1);
00384 
00385        c1 = cpositron[iZ][iT+1];
00386        c2 = cpositron[iZ+1][iT+1];
00387        cc2 = c1+ratZ*(c2-c1);
00388 
00389        corr = cc1+ratb2*(cc2-cc1);
00390 
00391        sigma *= sigmafactor/corr;
00392     }
00393   }
00394   else
00395   {
00396     c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
00397     c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
00398     if((AtomicNumber >= ZZ1) && (AtomicNumber <= ZZ2))
00399       sigma = c1+ratZ*(c2-c1) ;
00400     else if(AtomicNumber < ZZ1)
00401       sigma = AtomicNumber*AtomicNumber*c1/(ZZ1*ZZ1);
00402     else if(AtomicNumber > ZZ2)
00403       sigma = AtomicNumber*AtomicNumber*c2/(ZZ2*ZZ2);
00404   }
00405   return sigma;
00406 
00407 }
00408 
00409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00410 
00411 void G4UrbanMscModel96::StartTracking(G4Track* track)
00412 {
00413   SetParticle(track->GetDynamicParticle()->GetDefinition());
00414   firstStep = true; 
00415   inside = false;
00416   insideskin = false;
00417   tlimit = geombig;
00418   stepmin = tlimitminfix ;
00419   tlimitmin = 10.*stepmin ;
00420 }
00421 
00422 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00423 
00424 G4double G4UrbanMscModel96::ComputeTruePathLengthLimit(
00425                              const G4Track& track,
00426                              G4double& currentMinimalStep)
00427 {
00428   tPathLength = currentMinimalStep;
00429   const G4DynamicParticle* dp = track.GetDynamicParticle();
00430   
00431   G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
00432   G4StepStatus stepStatus = sp->GetStepStatus();
00433   couple = track.GetMaterialCutsCouple();
00434   SetCurrentCouple(couple); 
00435   currentMaterialIndex = couple->GetIndex();
00436   currentKinEnergy = dp->GetKineticEnergy();
00437   currentRange = GetRange(particle,currentKinEnergy,couple);
00438   lambda0 = GetTransportMeanFreePath(particle,currentKinEnergy);
00439   if(tPathLength > currentRange) { tPathLength = currentRange; }
00440 
00441   // stop here if small range particle
00442   if(inside || tPathLength < tlimitminfix) { 
00443     return ConvertTrueToGeom(tPathLength, currentMinimalStep); 
00444   }
00445   
00446   presafety = sp->GetSafety();
00447   /*
00448   G4cout << "G4Urban96::StepLimit tPathLength= " 
00449          <<tPathLength<<" safety= " << presafety
00450           << " range= " <<currentRange<< " lambda= "<<lambda0
00451          << " Alg: " << steppingAlgorithm <<G4endl;
00452   */
00453   // far from geometry boundary
00454   if(currentRange < presafety)
00455     {
00456       inside = true;
00457       return ConvertTrueToGeom(tPathLength, currentMinimalStep);  
00458     }
00459 
00460   // standard  version
00461   //
00462   if (steppingAlgorithm == fUseDistanceToBoundary)
00463     {
00464       //compute geomlimit and presafety 
00465       geomlimit = ComputeGeomLimit(track, presafety, currentRange);
00466 
00467       // is it far from boundary ?
00468       if(currentRange < presafety)
00469         {
00470           inside = true;
00471           return ConvertTrueToGeom(tPathLength, currentMinimalStep);   
00472         }
00473 
00474       smallstep += 1.;
00475       insideskin = false;
00476 
00477       if(firstStep || (stepStatus == fGeomBoundary))
00478         {
00479           rangeinit = currentRange;
00480           if(firstStep) smallstep = 1.e10;
00481           else  smallstep = 1.;
00482 
00483           //define stepmin here (it depends on lambda!)
00484           //rough estimation of lambda_elastic/lambda_transport
00485           G4double rat = currentKinEnergy/MeV ;
00486           rat = 1.e-3/(rat*(10.+rat)) ;
00487           //stepmin ~ lambda_elastic
00488           stepmin = rat*lambda0;
00489           skindepth = skin*stepmin;
00490           //define tlimitmin
00491           tlimitmin = 10.*stepmin;
00492           if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
00493           //G4cout << "rangeinit= " << rangeinit << " stepmin= " << stepmin
00494           //     << " tlimitmin= " << tlimitmin << " geomlimit= " << geomlimit <<G4endl;
00495           // constraint from the geometry
00496           if((geomlimit < geombig) && (geomlimit > geommin))
00497             {
00498               // geomlimit is a geometrical step length
00499               // transform it to true path length (estimation)
00500               if((1.-geomlimit/lambda0) > 0.)
00501                 geomlimit = -lambda0*log(1.-geomlimit/lambda0)+tlimitmin ;
00502 
00503               if(stepStatus == fGeomBoundary)
00504                 tgeom = geomlimit/facgeom;
00505               else
00506                 tgeom = 2.*geomlimit/facgeom;
00507             }
00508             else
00509               tgeom = geombig;
00510         }
00511 
00512 
00513       //step limit 
00514       tlimit = facrange*rangeinit;              
00515 
00516       //lower limit for tlimit
00517       if(tlimit < tlimitmin) tlimit = tlimitmin;
00518 
00519       if(tlimit > tgeom) tlimit = tgeom;
00520 
00521       //G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit  
00522       //      << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
00523 
00524       // shortcut
00525       if((tPathLength < tlimit) && (tPathLength < presafety) &&
00526          (smallstep > skin) && (tPathLength < geomlimit-0.999*skindepth))
00527         return ConvertTrueToGeom(tPathLength, currentMinimalStep);   
00528 
00529       // step reduction near to boundary
00530       if(smallstep <= skin)
00531         {
00532           tlimit = stepmin;
00533           insideskin = true;
00534         }
00535       else if(geomlimit < geombig)
00536         {
00537           if(geomlimit > skindepth)
00538             {
00539               if(tlimit > geomlimit-0.999*skindepth)
00540                 tlimit = geomlimit-0.999*skindepth;
00541             }
00542           else
00543             {
00544               insideskin = true;
00545               if(tlimit > stepmin) tlimit = stepmin;
00546             }
00547         }
00548 
00549       if(tlimit < stepmin) tlimit = stepmin;
00550 
00551       // randomize 1st step or 1st 'normal' step in volume
00552       if(firstStep || ((smallstep == skin) && !insideskin)) 
00553         { 
00554           G4double temptlimit = tlimit;
00555           if(temptlimit > tlimitmin)
00556           {
00557             do {
00558               temptlimit = G4RandGauss::shoot(tlimit,0.3*tlimit);        
00559                } while ((temptlimit < tlimitmin) || 
00560                         (temptlimit > 2.*tlimit-tlimitmin));
00561           }
00562           else
00563             temptlimit = tlimitmin;
00564           if(tPathLength > temptlimit) tPathLength = temptlimit;
00565         }
00566       else
00567         {  
00568           if(tPathLength > tlimit) tPathLength = tlimit  ; 
00569         }
00570 
00571     }
00572     // for 'normal' simulation with or without magnetic field 
00573     //  there no small step/single scattering at boundaries
00574   else if(steppingAlgorithm == fUseSafety)
00575     {
00576       // compute presafety again if presafety <= 0 and no boundary
00577       // i.e. when it is needed for optimization purposes
00578       if((stepStatus != fGeomBoundary) && (presafety < tlimitminfix)) 
00579         presafety = ComputeSafety(sp->GetPosition(),tPathLength); 
00580       /*
00581       G4cout << "presafety= " << presafety
00582              << " firstStep= " << firstStep
00583              << " stepStatus= " << stepStatus 
00584              << G4endl;
00585       */
00586       // is far from boundary
00587       if(currentRange < presafety)
00588         {
00589           inside = true;
00590           return ConvertTrueToGeom(tPathLength, currentMinimalStep);  
00591         }
00592 
00593       if(firstStep || stepStatus == fGeomBoundary)
00594       {
00595         rangeinit = currentRange;
00596         fr = facrange;
00597         // 9.1 like stepping for e+/e- only (not for muons,hadrons)
00598         if(mass < masslimite) 
00599         {
00600           if(lambda0 > currentRange)
00601             rangeinit = lambda0;
00602           if(lambda0 > lambdalimit)
00603             fr *= 0.75+0.25*lambda0/lambdalimit;
00604         }
00605 
00606         //lower limit for tlimit
00607         G4double rat = currentKinEnergy/MeV ;
00608         rat = 1.e-3/(rat*(10.+rat)) ;
00609         stepmin = lambda0*rat;
00610         tlimitmin = 10.*stepmin;
00611         if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
00612       }
00613       //step limit
00614       tlimit = fr*rangeinit;               
00615 
00616       if(tlimit < facsafety*presafety)
00617         tlimit = facsafety*presafety;
00618 
00619       //lower limit for tlimit
00620       if(tlimit < tlimitmin) tlimit = tlimitmin;
00621       
00622       if(firstStep || stepStatus == fGeomBoundary)
00623       {
00624         G4double temptlimit = tlimit;
00625         if(temptlimit > tlimitmin)
00626         {
00627           do {
00628             temptlimit = G4RandGauss::shoot(tlimit,0.3*tlimit);
00629              } while ((temptlimit < tlimitmin) ||
00630                       (temptlimit > 2.*tlimit-tlimitmin));
00631         }
00632         else
00633           temptlimit = tlimitmin;
00634 
00635         if(tPathLength > temptlimit) tPathLength = temptlimit;
00636       }
00637       else
00638         if(tPathLength > tlimit) tPathLength = tlimit;
00639     }
00640   
00641   // version similar to 7.1 (needed for some experiments)
00642   else
00643     {
00644       if (stepStatus == fGeomBoundary)
00645         {
00646           if (currentRange > lambda0) tlimit = facrange*currentRange;
00647           else                        tlimit = facrange*lambda0;
00648 
00649           if(tlimit < tlimitmin) tlimit = tlimitmin;
00650           if(tPathLength > tlimit) tPathLength = tlimit;
00651         }
00652     }
00653   //G4cout << "tPathLength= " << tPathLength 
00654   //     << " currentMinimalStep= " << currentMinimalStep << G4endl;
00655   return ConvertTrueToGeom(tPathLength, currentMinimalStep);
00656 }
00657 
00658 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00659 
00660 G4double G4UrbanMscModel96::ComputeGeomPathLength(G4double)
00661 {
00662   firstStep = false; 
00663   lambdaeff = lambda0;
00664   par1 = -1. ;  
00665   par2 = par3 = 0. ;  
00666 
00667   //  do the true -> geom transformation
00668   zPathLength = tPathLength;
00669 
00670   // z = t for very small tPathLength
00671   if(tPathLength < tlimitminfix) return zPathLength;
00672 
00673   // this correction needed to run MSC with eIoni and eBrem inactivated
00674   // and makes no harm for a normal run
00675   // It is already checked
00676   // if(tPathLength > currentRange)
00677   //  tPathLength = currentRange ;
00678 
00679   G4double tau   = tPathLength/lambda0 ;
00680 
00681   if ((tau <= tausmall) || insideskin) {
00682     zPathLength  = tPathLength;
00683     if(zPathLength > lambda0) zPathLength = lambda0;
00684     return zPathLength;
00685   }
00686 
00687   G4double zmean = tPathLength;
00688   if (tPathLength < currentRange*dtrl) {
00689     if(tau < taulim) zmean = tPathLength*(1.-0.5*tau) ;
00690     else             zmean = lambda0*(1.-exp(-tau));
00691     zPathLength = zmean ;
00692     return zPathLength;    
00693 
00694   } else if(currentKinEnergy < mass || tPathLength == currentRange)  {
00695     par1 = 1./currentRange ;
00696     par2 = 1./(par1*lambda0) ;
00697     par3 = 1.+par2 ;
00698     if(tPathLength < currentRange)
00699       zmean = (1.-exp(par3*log(1.-tPathLength/currentRange)))/(par1*par3) ;
00700     else {
00701       zmean = 1./(par1*par3) ;
00702     }
00703     zPathLength = zmean ;
00704     return zPathLength;    
00705 
00706   } else {
00707     G4double T1 = GetEnergy(particle,currentRange-tPathLength,couple);
00708     G4double lambda1 = GetTransportMeanFreePath(particle,T1);
00709 
00710     par1 = (lambda0-lambda1)/(lambda0*tPathLength);
00711     par2 = 1./(par1*lambda0);
00712     par3 = 1.+par2 ;
00713     zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3);
00714   }
00715 
00716   zPathLength = zmean;
00717 
00718   //  sample z
00719   if(samplez)
00720   {
00721     const G4double  ztmax = 0.999 ;
00722     G4double zt = zmean/tPathLength ;
00723 
00724     if (tPathLength > stepmin && zt < ztmax)              
00725     {
00726       G4double u,cz1;
00727       if(zt >= third)
00728       {
00729         G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
00730         cz1 = 1.+cz ;
00731         G4double u0 = cz/cz1 ;
00732         G4double grej ;
00733         do {
00734             u = exp(log(G4UniformRand())/cz1) ;
00735             grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
00736            } while (grej < G4UniformRand()) ;
00737       }
00738       else
00739       {
00740         u = 2.*zt*G4UniformRand();
00741       }
00742       zPathLength = tPathLength*u ;
00743     }
00744   }
00745 
00746   if(zPathLength > lambda0) { zPathLength = lambda0; }
00747   //G4cout<< "zPathLength= "<< zPathLength<< " lambda1= " << lambda0 << G4endl;
00748   return zPathLength;
00749 }
00750 
00751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00752 
00753 G4double G4UrbanMscModel96::ComputeTrueStepLength(G4double geomStepLength)
00754 {
00755   // step defined other than transportation 
00756   if(geomStepLength == zPathLength)
00757     { return tPathLength; }
00758 
00759   zPathLength = geomStepLength;
00760 
00761   // t = z for very small step
00762   if(geomStepLength < tlimitminfix) { 
00763     tPathLength = geomStepLength; 
00764   
00765   // recalculation
00766   } else {
00767 
00768     G4double tlength = geomStepLength;
00769     if((geomStepLength > lambda0*tausmall) && !insideskin) {
00770 
00771       if(par1 <  0.) {
00772         tlength = -lambda0*log(1.-geomStepLength/lambda0) ;
00773       } else {
00774         if(par1*par3*geomStepLength < 1.) {
00775           tlength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
00776         } else {
00777           tlength = currentRange;
00778         }
00779       }
00780       if(tlength < geomStepLength)   { tlength = geomStepLength; }
00781       else if(tlength > tPathLength) { tlength = tPathLength; }
00782     }  
00783     tPathLength = tlength; 
00784   }
00785   //G4cout << "Urban96::ComputeTrueLength: tPathLength= " << tPathLength 
00786   //     << " step= " << geomStepLength << G4endl;
00787 
00788   return tPathLength;
00789 }
00790 
00791 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00792 
00793 G4double G4UrbanMscModel96::ComputeTheta0(G4double trueStepLength,
00794                                           G4double KineticEnergy)
00795 {
00796   // for all particles take the width of the central part
00797   //  from a  parametrization similar to the Highland formula
00798   // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
00799   static const G4double c_highland = 13.6*MeV ;
00800   G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
00801                          KineticEnergy*(KineticEnergy+2.*mass)/
00802                       ((currentKinEnergy+mass)*(KineticEnergy+mass)));
00803   y = trueStepLength/currentRadLength;
00804   G4double theta0 = c_highland*std::abs(charge)*sqrt(y)/betacp;
00805   y = log(y);
00806   // correction factor from e- scattering data
00807   G4double corr = coeffth1+coeffth2*y;                
00808 
00809   theta0 *= corr ;                                               
00810 
00811   return theta0;
00812 }
00813 
00814 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00815 
00816 G4ThreeVector& 
00817 G4UrbanMscModel96::SampleScattering(const G4ThreeVector& oldDirection,
00818                                     G4double safety)
00819 {
00820   fDisplacement.set(0.0,0.0,0.0);
00821   G4double kineticEnergy = currentKinEnergy;
00822   if (tPathLength > currentRange*dtrl) {
00823     kineticEnergy = GetEnergy(particle,currentRange-tPathLength,couple);
00824   } else {
00825     kineticEnergy -= tPathLength*GetDEDX(particle,currentKinEnergy,couple);
00826   }
00827 
00828   if((kineticEnergy <= eV) || (tPathLength <= tlimitminfix) ||
00829      (tPathLength/tausmall < lambda0)) { return fDisplacement; }
00830 
00831   G4double cth  = SampleCosineTheta(tPathLength,kineticEnergy);
00832 
00833   // protection against 'bad' cth values
00834   if(std::fabs(cth) > 1.) { return fDisplacement; }
00835 
00836   // extra protection agaist high energy particles backscattered 
00837     //G4cout << "Warning: large scattering E(MeV)= " << kineticEnergy 
00838     //     << " s(mm)= " << tPathLength/mm
00839     //     << " 1-cosTheta= " << 1.0 - cth << G4endl;
00840     // do Gaussian central scattering
00841   //  if(kineticEnergy > 5*GeV && cth < 0.9) {
00842   /*
00843   if(cth < 1.0 - 1000*tPathLength/lambda0 
00844      && cth < 0.8 && kineticEnergy > 20*MeV) { 
00845     G4ExceptionDescription ed;
00846     ed << particle->GetParticleName()
00847        << " E(MeV)= " << kineticEnergy/MeV
00848        << " Step(mm)= " << tPathLength/mm
00849        << " in " << CurrentCouple()->GetMaterial()->GetName()
00850        << " CosTheta= " << cth 
00851        << " is too big";
00852     G4Exception("G4UrbanMscModel96::SampleScattering","em0004",
00853                 JustWarning, ed,"");
00854 
00855     if(kineticEnergy > GeV && cth < 0.0) {
00856     do {
00857       cth = 1.0 + 2*log(G4UniformRand())*tPathLength/lambda0;
00858     } while(cth < -1.0);
00859   }
00860 */
00861 
00862   G4double sth  = sqrt((1.0 - cth)*(1.0 + cth));
00863   G4double phi  = twopi*G4UniformRand();
00864   G4double dirx = sth*cos(phi);
00865   G4double diry = sth*sin(phi);
00866 
00867   G4ThreeVector newDirection(dirx,diry,cth);
00868   newDirection.rotateUz(oldDirection);
00869   fParticleChange->ProposeMomentumDirection(newDirection);
00870   /*
00871   G4cout << "G4UrbanMscModel96::SampleSecondaries: e(MeV)= " << kineticEnergy
00872          << " sinTheta= " << sth << " safety(mm)= " << safety
00873          << " trueStep(mm)= " << tPathLength
00874          << " geomStep(mm)= " << zPathLength
00875          << G4endl;
00876   */
00877   if (latDisplasment && safety > tlimitminfix) {
00878 
00879     G4double r = SampleDisplacement();
00880     /*    
00881     G4cout << "G4UrbanMscModel96::SampleSecondaries: e(MeV)= " << kineticEnergy
00882            << " sinTheta= " << sth << " r(mm)= " << r
00883            << " trueStep(mm)= " << tPathLength
00884            << " geomStep(mm)= " << zPathLength
00885            << G4endl;
00886     */
00887     if(r > 0.)
00888       {
00889         G4double latcorr = LatCorrelation();
00890         if(latcorr > r) latcorr = r;
00891 
00892         // sample direction of lateral displacement
00893         // compute it from the lateral correlation
00894         G4double Phi = 0.;
00895         if(std::abs(r*sth) < latcorr)
00896           Phi  = twopi*G4UniformRand();
00897         else
00898         {
00899           G4double psi = std::acos(latcorr/(r*sth));
00900           if(G4UniformRand() < 0.5)
00901             Phi = phi+psi;
00902           else
00903             Phi = phi-psi;
00904         }
00905 
00906         dirx = std::cos(Phi);
00907         diry = std::sin(Phi);
00908 
00909         fDisplacement.set(r*dirx,r*diry,0.0);
00910         fDisplacement.rotateUz(oldDirection);
00911       }
00912   }
00913   return fDisplacement;
00914 }
00915 
00916 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00917 
00918 G4double G4UrbanMscModel96::SampleCosineTheta(G4double trueStepLength,
00919                                               G4double KineticEnergy)
00920 {
00921   G4double cth = 1. ;
00922   G4double tau = trueStepLength/lambda0;
00923   currentTau   = tau;
00924   lambdaeff    = lambda0;
00925 
00926   Zeff = couple->GetMaterial()->GetTotNbOfElectPerVolume()/
00927          couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
00928 
00929   if(Zold != Zeff)  
00930     UpdateCache();
00931 
00932   G4double lambda1 = GetTransportMeanFreePath(particle,KineticEnergy);
00933   if(std::fabs(lambda1/lambda0 - 1) > 0.01 && lambda1 > 0.)
00934   {
00935     // mean tau value
00936     tau = trueStepLength*log(lambda0/lambda1)/(lambda0-lambda1);
00937   }
00938 
00939   currentTau = tau ;
00940   lambdaeff = trueStepLength/currentTau;
00941   currentRadLength = couple->GetMaterial()->GetRadlen();
00942 
00943   if (tau >= taubig) { cth = -1.+2.*G4UniformRand(); }
00944   else if (tau >= tausmall) {
00945     static const G4double numlim = 0.01;
00946     G4double xmeanth, x2meanth;
00947     if(tau < numlim) {
00948       xmeanth = 1.0 - tau*(1.0 - 0.5*tau);
00949       x2meanth= 1.0 - tau*(5.0 - 6.25*tau)/3.;
00950     } else {
00951       xmeanth = exp(-tau);
00952       x2meanth = (1.+2.*exp(-2.5*tau))/3.;
00953     }
00954     G4double relloss = 1.-KineticEnergy/currentKinEnergy;
00955 
00956     if(relloss > rellossmax) 
00957       return SimpleScattering(xmeanth,x2meanth);
00958 
00959     // is step extreme small ?
00960     G4bool extremesmallstep = false ;
00961     G4double tsmall = tlimitmin ;
00962     G4double theta0 = 0.;
00963     if(trueStepLength > tsmall) {
00964        theta0 = ComputeTheta0(trueStepLength,KineticEnergy);
00965     } else {
00966       G4double rate = trueStepLength/tsmall ;
00967       if(G4UniformRand() < rate) {
00968         theta0 = ComputeTheta0(tsmall,KineticEnergy);
00969         extremesmallstep = true ;
00970       }
00971     }
00972     //G4cout << "Theta0= " << theta0 << " theta0max= " << theta0max 
00973     //       << "  sqrt(tausmall)= " << sqrt(tausmall) << G4endl;
00974 
00975     // protection for very small angles
00976     G4double theta2 = theta0*theta0;
00977 
00978     if(theta2 < tausmall) { return cth; }
00979     
00980     if(theta0 > theta0max) {
00981       return SimpleScattering(xmeanth,x2meanth);
00982     }
00983 
00984     G4double x = theta2*(1.0 - theta2/12.);
00985     if(theta2 > numlim) {
00986       G4double sth = 2*sin(0.5*theta0);
00987       x = sth*sth;
00988     }
00989 
00990     // parameter for tail
00991     G4double ltau= log(tau);
00992     G4double u   = exp(ltau/6.);
00993     if(extremesmallstep)  u = exp(log(tsmall/lambda0)/6.);
00994     G4double xx  = log(lambdaeff/currentRadLength);
00995     G4double xsi = coeffc1+u*(coeffc2+coeffc3*u)+coeffc4*xx;
00996 
00997     // tail should not be too big
00998     if(xsi < 1.9) { 
00999       /*
01000       if(KineticEnergy > 20*MeV && xsi < 1.6) {
01001         G4cout << "G4UrbanMscModel96::SampleCosineTheta: E(GeV)= " 
01002                << KineticEnergy/GeV 
01003                << " !!** c= " << xsi
01004                << " **!! length(mm)= " << trueStepLength << " Zeff= " << Zeff 
01005                << " " << couple->GetMaterial()->GetName()
01006                << " tau= " << tau << G4endl;
01007       }
01008       */
01009       xsi = 1.9; 
01010     }
01011 
01012     G4double c = xsi;
01013 
01014     if(fabs(c-3.) < 0.001)      { c = 3.001; }
01015     else if(fabs(c-2.) < 0.001) { c = 2.001; }
01016 
01017     G4double c1 = c-1.;
01018 
01019     G4double ea = exp(-xsi);
01020     G4double eaa = 1.-ea ;
01021     G4double xmean1 = 1.-(1.-(1.+xsi)*ea)*x/eaa;
01022     G4double x0 = 1. - xsi*x;
01023 
01024     // G4cout << " xmean1= " << xmean1 << "  xmeanth= " << xmeanth << G4endl;
01025 
01026     if(xmean1 <= 0.999*xmeanth) {
01027       return SimpleScattering(xmeanth,x2meanth);
01028     }
01029     //from continuity of derivatives
01030     G4double b = 1.+(c-xsi)*x;
01031 
01032     G4double b1 = b+1.;
01033     G4double bx = c*x;
01034 
01035     G4double eb1 = pow(b1,c1);
01036     G4double ebx = pow(bx,c1);
01037     G4double d = ebx/eb1;
01038 
01039     G4double xmean2 = (x0 + d - (bx - b1*d)/(c-2.))/(1. - d);
01040 
01041     G4double f1x0 = ea/eaa;
01042     G4double f2x0 = c1/(c*(1. - d));
01043     G4double prob = f2x0/(f1x0+f2x0);
01044 
01045     G4double qprob = xmeanth/(prob*xmean1+(1.-prob)*xmean2);
01046 
01047     // sampling of costheta
01048     //G4cout << "c= " << c << " qprob= " << qprob << " eb1= " << eb1
01049     // << " c1= " << c1 << " b1= " << b1 << " bx= " << bx << " eb1= " << eb1
01050     //       << G4endl;
01051     if(G4UniformRand() < qprob)
01052     {
01053       G4double var = 0;
01054       if(G4UniformRand() < prob) {
01055         cth = 1.+log(ea+G4UniformRand()*eaa)*x;
01056       } else {
01057         var = (1.0 - d)*G4UniformRand();
01058         if(var < numlim*d) {
01059           var /= (d*c1); 
01060           cth = -1.0 + var*(1.0 - 0.5*var*c)*(2. + (c - xsi)*x);
01061         } else {
01062           cth = 1. + x*(c - xsi - c*pow(var + d, -1.0/c1));
01063           //b-b1*bx/exp(log(ebx+(eb1-ebx)*G4UniformRand())/c1) ;
01064         }
01065       } 
01066       if(KineticEnergy > 5*GeV && cth < 0.9) {
01067         G4cout << "G4UrbanMscModel96::SampleCosineTheta: E(GeV)= " 
01068                << KineticEnergy/GeV 
01069                << " 1-cosT= " << 1 - cth
01070                << " length(mm)= " << trueStepLength << " Zeff= " << Zeff 
01071                << " tau= " << tau
01072                << " prob= " << prob << " var= " << var << G4endl;
01073         G4cout << "  c= " << c << " qprob= " << qprob << " eb1= " << eb1
01074                << " ebx= " << ebx
01075                << " c1= " << c1 << " b= " << b << " b1= " << b1 
01076                << " bx= " << bx << " d= " << d
01077                << " ea= " << ea << " eaa= " << eaa << G4endl;
01078       }
01079     }
01080     else {
01081       cth = -1.+2.*G4UniformRand();
01082       if(KineticEnergy > 5*GeV) {
01083         G4cout << "G4UrbanMscModel96::SampleCosineTheta: E(GeV)= " 
01084                << KineticEnergy/GeV 
01085                << " length(mm)= " << trueStepLength << " Zeff= " << Zeff 
01086                << " qprob= " << qprob << G4endl;
01087       }
01088     }
01089   }
01090   return cth ;
01091 }
01092 
01093 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
01094 
01095 G4double G4UrbanMscModel96::SimpleScattering(G4double xmeanth,G4double x2meanth)
01096 {
01097   // 'large angle scattering'
01098   // 2 model functions with correct xmean and x2mean
01099   G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
01100   G4double prob = (a+2.)*xmeanth/a;
01101 
01102   // sampling
01103   G4double cth = 1.;
01104   if(G4UniformRand() < prob)
01105     cth = -1.+2.*exp(log(G4UniformRand())/(a+1.));
01106   else
01107     cth = -1.+2.*G4UniformRand();
01108   return cth;
01109 }
01110 
01111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
01112 
01113 G4double G4UrbanMscModel96::SampleDisplacement()
01114 {
01115   G4double r = 0.0;
01116   if ((currentTau >= tausmall) && !insideskin) {
01117     G4double rmax = sqrt((tPathLength-zPathLength)*(tPathLength+zPathLength));
01118     r = rmax*exp(log(G4UniformRand())/3.);
01119   }
01120   return r;
01121 }
01122 
01123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
01124 
01125 G4double G4UrbanMscModel96::LatCorrelation()
01126 {
01127   static const G4double kappa = 2.5;
01128   static const G4double kappami1 = kappa-1.;
01129 
01130   G4double latcorr = 0.;
01131   if((currentTau >= tausmall) && !insideskin)
01132   {
01133     if(currentTau < taulim)
01134       latcorr = lambdaeff*kappa*currentTau*currentTau*
01135                 (1.-(kappa+1.)*currentTau/3.)/3.;
01136     else
01137     {
01138       G4double etau = 0.;
01139       if(currentTau < taubig) etau = exp(-currentTau);
01140       latcorr = -kappa*currentTau;
01141       latcorr = exp(latcorr)/kappami1;
01142       latcorr += 1.-kappa*etau/kappami1 ;
01143       latcorr *= 2.*lambdaeff/3. ;
01144     }
01145   }
01146 
01147   return latcorr;
01148 }
01149 
01150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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