processes/electromagnetic/standard/src/G4MscModel71.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: G4MscModel71.cc,v 1.5 2007/05/22 17:34:36 vnivanch Exp $
00027 // GEANT4 tag $Name: geant4-09-01-patch-01 $
00028 //
00029 // -------------------------------------------------------------------
00030 //
00031 // GEANT4 Class file
00032 //
00033 //
00034 // File name:   G4MscModel71
00035 //
00036 // Author:      Laszlo Urban
00037 //
00038 // Creation date: 03.03.2001
00039 //
00040 // Modifications:
00041 //
00042 // 27-03-03 Move model part from G4MultipleScattering (V.Ivanchenko)
00043 // 23-05-03 important change in angle distribution for muons/hadrons
00044 //          the central part now is similar to the Highland parametrization +
00045 //          minor correction in angle sampling algorithm (for all particles)
00046 //          (L.Urban)
00047 // 30-05-03 misprint in SampleCosineTheta corrected(L.Urban)
00048 // 27-03-03 Rename (V.Ivanchenko)
00049 // 05-08-03 angle distribution has been modified (L.Urban)
00050 // 06-11-03 precision problems solved for high energy (PeV) particles
00051 //          change in the tail of the angular distribution
00052 //          highKinEnergy is set to 100 PeV (L.Urban) 
00053 //
00054 // 10-11-03 highKinEnergy is set back to 100 TeV, some tail tuning +
00055 //          cleaning (L.Urban) 
00056 // 26-11-03 correction in TrueStepLength : 
00057 //          trueLength <= currentRange (L.Urban) 
00058 // 01-03-04 signature changed in SampleCosineTheta,
00059 //          energy dependence calculations has been simplified,
00060 // 11-03-04 corrections in GeomPathLength,TrueStepLength,
00061 //          SampleCosineTheta
00062 // 23-04-04 true -> geom and geom -> true transformation has been
00063 //          rewritten, changes in the angular distribution (L.Urban)
00064 // 19-07-04 correction in SampleCosineTheta in order to avoid
00065 //          num. precision problems at high energy/small step(L.Urban) 
00066 // 17-08-04 changes in the angle distribution (slightly modified
00067 //          Highland formula for the width of the central part,
00068 //          changes in the numerical values of some other parameters)
00069 //          ---> approximately step independent distribution (L.Urban)
00070 // 21-09-04 change in the tail of the angular distribution (L.Urban)
00071 //
00072 // 03-11-04 precision problem for very high energy ions and small stepsize
00073 //          solved in SampleCosineTheta (L.Urban).
00074 // 15-04-05 optimize internal interface - add SampleSecondaries method (V.Ivanchenko)
00075 // 03-10-05 Model is freezed with the name McsModel71 (V.Ivanchenko)
00076 // 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
00077 //
00078 
00079 // Class Description:
00080 //
00081 // Implementation of the model of multiple scattering based on
00082 // H.W.Lewis Phys Rev 78 (1950) 526 and others
00083 
00084 // -------------------------------------------------------------------
00085 //
00086 
00087 
00088 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00089 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00090 
00091 #include "G4MscModel71.hh"
00092 #include "Randomize.hh"
00093 #include "G4Electron.hh"
00094 #include "G4LossTableManager.hh"
00095 #include "G4PhysicsTable.hh"
00096 #include "G4ParticleChangeForMSC.hh"
00097 #include "G4TransportationManager.hh"
00098 #include "G4Navigator.hh"
00099 
00100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00101 
00102 using namespace std;
00103 
00104 G4MscModel71::G4MscModel71(G4double& m_dtrl, G4double& m_NuclCorrPar,
00105                            G4double& m_FactPar, G4double& m_factail,
00106                        G4bool& m_samplez, const G4String& nam)
00107   : G4VEmModel(nam),
00108   taubig(8.0),
00109   tausmall(1.e-20),
00110   taulim(1.e-6),
00111   dtrl(m_dtrl),
00112   NuclCorrPar (m_NuclCorrPar),
00113   FactPar(m_FactPar),
00114   factail(m_factail),
00115   samplez(m_samplez),
00116   isInitialized(false)
00117 {
00118   stepmin       = 1.e-6*mm;
00119   currentRange  = 0. ;
00120 }
00121 
00122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00123 
00124 G4MscModel71::~G4MscModel71()
00125 {}
00126 
00127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00128 
00129 void G4MscModel71::Initialise(const G4ParticleDefinition* p,
00130                               const G4DataVector&)
00131 {
00132   if(isInitialized) return;
00133   // set values of some data members
00134   sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
00135   particle = p;
00136   mass = particle->GetPDGMass();
00137   charge = particle->GetPDGCharge()/eplus;
00138   b = 1. ;
00139   xsi = 3.00 ;
00140 
00141   if(pParticleChange)
00142     fParticleChange = reinterpret_cast<G4ParticleChangeForMSC*>(pParticleChange);
00143   else
00144     fParticleChange = new G4ParticleChangeForMSC();
00145 
00146   navigator = G4TransportationManager::GetTransportationManager()
00147     ->GetNavigatorForTracking();
00148 }
00149 
00150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00151 
00152 G4double G4MscModel71::ComputeCrossSectionPerAtom(
00153                              const G4ParticleDefinition* part,
00154                                    G4double KineticEnergy,
00155                                    G4double AtomicNumber,
00156                                    G4double AtomicWeight, 
00157                                    G4double,
00158                                    G4double)
00159 {
00160   const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
00161                              Bohr_radius*Bohr_radius/(hbarc*hbarc);
00162   const G4double epsmin = 1.e-4 , epsmax = 1.e10;
00163 
00164   const G4double Zdat[15] = { 4., 6.,13.,20.,26.,29.,32.,38.,47.,
00165                              50.,56.,64.,74.,79.,82. };
00166 
00167   const G4double Tdat[23] = {0.0001*MeV,0.0002*MeV,0.0004*MeV,0.0007*MeV,
00168                              0.001*MeV,0.002*MeV,0.004*MeV,0.007*MeV,
00169                              0.01*MeV,0.02*MeV,0.04*MeV,0.07*MeV,
00170                              0.1*MeV,0.2*MeV,0.4*MeV,0.7*MeV,
00171                              1.*MeV,2.*MeV,4.*MeV,7.*MeV,10.*MeV,20.*MeV,
00172                              10000.0*MeV};
00173 
00174  // corr. factors for e-/e+ lambda
00175           G4double celectron[15][23] =
00176           {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
00177             1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
00178             1.112,1.108,1.100,1.093,1.089,1.087,0.7235     },
00179            {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
00180             1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
00181             1.109,1.105,1.097,1.090,1.086,1.082,0.7925     },
00182            {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
00183             1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
00184             1.131,1.124,1.113,1.104,1.099,1.098,0.9147     },
00185            {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
00186             1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
00187             1.112,1.105,1.096,1.089,1.085,1.098,0.9700     },
00188            {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
00189             1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
00190             1.073,1.070,1.064,1.059,1.056,1.056,1.0022     },
00191            {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
00192             1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
00193             1.074,1.070,1.063,1.059,1.056,1.052,1.0158     },
00194            {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
00195             1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
00196             1.068,1.064,1.059,1.054,1.051,1.050,1.0284     },
00197            {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
00198             1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
00199             1.039,1.037,1.034,1.031,1.030,1.036,1.0515     },
00200            {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
00201             1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
00202             1.031,1.028,1.024,1.022,1.021,1.024,1.0834     },
00203            {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
00204             1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
00205             1.020,1.017,1.015,1.013,1.013,1.020,1.0937     },
00206            {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
00207             1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
00208             0.995,0.993,0.993,0.993,0.993,1.011,1.1140     },
00209            {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
00210             1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
00211             0.974,0.972,0.973,0.974,0.975,0.987,1.1410     },
00212            {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
00213             1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
00214             0.950,0.947,0.949,0.952,0.954,0.963,1.1750     },
00215            {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
00216             1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
00217             0.941,0.938,0.940,0.944,0.946,0.954,1.1922     },
00218            {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
00219             1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
00220             0.933,0.930,0.933,0.936,0.939,0.949,1.2026     }};
00221            G4double cpositron[15][23] = {
00222            {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
00223             1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
00224             1.131,1.126,1.117,1.108,1.103,1.100,0.7235     },
00225            {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
00226             1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
00227             1.138,1.132,1.122,1.113,1.108,1.102,0.7925     },
00228            {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
00229             1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
00230             1.203,1.190,1.173,1.159,1.151,1.145,0.9147     },
00231            {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
00232             1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
00233             1.225,1.210,1.191,1.175,1.166,1.174,0.9700     },
00234            {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
00235             1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
00236             1.217,1.203,1.184,1.169,1.160,1.151,1.0022     },
00237            {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
00238             1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
00239             1.237,1.222,1.201,1.184,1.174,1.159,1.0158     },
00240            {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
00241             1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
00242             1.252,1.234,1.212,1.194,1.183,1.170,1.0284     },
00243            {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
00244             2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
00245             1.254,1.237,1.214,1.195,1.185,1.179,1.0515     },
00246            {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
00247             2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
00248             1.312,1.288,1.258,1.235,1.221,1.205,1.0834     },
00249            {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
00250             2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
00251             1.320,1.294,1.264,1.240,1.226,1.214,1.0937     },
00252            {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
00253             2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
00254             1.328,1.302,1.270,1.245,1.231,1.233,1.1140     },
00255            {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
00256             2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
00257             1.361,1.330,1.294,1.267,1.251,1.239,1.1410     },
00258            {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
00259             3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
00260             1.409,1.372,1.330,1.298,1.280,1.258,1.1750     },
00261            {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
00262             3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
00263             1.442,1.400,1.354,1.319,1.299,1.272,1.1922     },
00264            {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
00265             3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
00266             1.456,1.412,1.364,1.328,1.307,1.282,1.2026     }};
00267 
00268   G4double sigma;
00269   if (part != particle ) {
00270     particle = part;
00271     mass = particle->GetPDGMass();
00272     charge = particle->GetPDGCharge()/eplus;
00273   }
00274 
00275   G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
00276 
00277   // correction if particle .ne. e-/e+
00278   // compute equivalent kinetic energy
00279   // lambda depends on p*beta ....
00280 
00281   G4double eKineticEnergy = KineticEnergy;
00282 
00283   if((particle->GetParticleName() != "e-") &&
00284      (particle->GetParticleName() != "e+") )
00285   {
00286      G4double TAU = KineticEnergy/mass ;
00287      G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
00288      G4double w = c-2. ;
00289      G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
00290      eKineticEnergy = electron_mass_c2*tau ;
00291   }
00292 
00293   G4double ChargeSquare = charge*charge;
00294 
00295   G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
00296   G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00297                                  /(eTotalEnergy*eTotalEnergy);
00298   G4double bg2   = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
00299                                  /(electron_mass_c2*electron_mass_c2);
00300 
00301   G4double eps = epsfactor*bg2/Z23;
00302 
00303   if     (eps<epsmin)  sigma = 2.*eps*eps;
00304   else if(eps<epsmax)  sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
00305   else                 sigma = log(2.*eps)-1.+1./eps;
00306 
00307   sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
00308 
00309   // nuclear size effect correction for high energy
00310   // ( a simple approximation at present)
00311   G4double corrnuclsize,a,w1,w2,w;
00312 
00313   G4double x0 = 1. - NuclCorrPar*mass/(KineticEnergy*
00314                exp(log(AtomicWeight/(g/mole))/3.));
00315   if ( x0 < -1. || eKineticEnergy  <= 10.*MeV)
00316       {
00317         x0 = -1.;
00318         corrnuclsize = 1.;
00319       }
00320   else
00321       {
00322         a = 1.+1./eps;
00323         if (eps > epsmax) w1=log(2.*eps)+1./eps-3./(8.*eps*eps);
00324         else              w1=log((a+1.)/(a-1.))-2./(a+1.);
00325         w = 1./((1.-x0)*eps);
00326         if (w < epsmin)   w2=-log(w)-1.+2.*w-1.5*w*w;
00327         else              w2 = log((a-x0)/(a-1.))-(1.-x0)/(a-x0);
00328         corrnuclsize = w1/w2;
00329         corrnuclsize = exp(-FactPar*mass/KineticEnergy)*
00330                       (corrnuclsize-1.)+1.;
00331       }
00332 
00333   // interpolate in AtomicNumber and beta2
00334   // get bin number in Z
00335   G4int iZ = 14;
00336   while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
00337   if (iZ==14)                               iZ = 13;
00338   if (iZ==-1)                               iZ = 0 ;
00339 
00340   G4double Z1 = Zdat[iZ];
00341   G4double Z2 = Zdat[iZ+1];
00342   G4double ratZ = (AtomicNumber-Z1)/(Z2-Z1);
00343 
00344   // get bin number in T (beta2)
00345   G4int iT = 22;
00346   while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
00347   if(iT==22)                                  iT = 21;
00348   if(iT==-1)                                  iT = 0 ;
00349 
00350   //  calculate betasquare values
00351   G4double T = Tdat[iT],   E = T + electron_mass_c2;
00352   G4double b2small = T*(E+electron_mass_c2)/(E*E);
00353   T = Tdat[iT+1]; E = T + electron_mass_c2;
00354   G4double b2big = T*(E+electron_mass_c2)/(E*E);
00355   G4double ratb2 = (beta2-b2small)/(b2big-b2small);
00356   G4double c1,c2,cc1,cc2,corr;
00357   if (charge < 0.)
00358     {
00359        c1 = celectron[iZ][iT];
00360        c2 = celectron[iZ+1][iT];
00361        cc1 = c1+ratZ*(c2-c1);
00362 
00363        c1 = celectron[iZ][iT+1];
00364        c2 = celectron[iZ+1][iT+1];
00365        cc2 = c1+ratZ*(c2-c1);
00366 
00367        corr = cc1+ratb2*(cc2-cc1);
00368        sigma /= corr;
00369     }
00370 
00371   if (charge > 0.)
00372     {
00373        c1 = cpositron[iZ][iT];
00374        c2 = cpositron[iZ+1][iT];
00375        cc1 = c1+ratZ*(c2-c1);
00376 
00377        c1 = cpositron[iZ][iT+1];
00378        c2 = cpositron[iZ+1][iT+1];
00379        cc2 = c1+ratZ*(c2-c1);
00380 
00381        corr = cc1+ratb2*(cc2-cc1);
00382        sigma /= corr;
00383     }
00384 
00385   sigma *= sigmafactor;
00386 
00387   //  nucl. size correction for particles other than e+/e- only at present !!!!
00388   if((particle->GetParticleName() != "e-") &&
00389      (particle->GetParticleName() != "e+")   )
00390      sigma /= corrnuclsize;
00391   //  G4cout << "e= " << KineticEnergy << " sigma= " << sigma << G4endl;
00392 
00393   return sigma;
00394 }
00395 
00396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00397 
00398 G4double G4MscModel71::GeomPathLength(
00399                           G4PhysicsTable* theLambdaTable,
00400                     const G4MaterialCutsCouple* couple,
00401                     const G4ParticleDefinition* theParticle,
00402                           G4double& T0,
00403                           G4double lambda,
00404                           G4double range,
00405                           G4double truePathLength)
00406 {
00407   //  do the true -> geom transformation
00408   const G4double  ztmax = 101./103. ;
00409   if (theParticle != particle ) {
00410     particle = theParticle;
00411     mass = particle->GetPDGMass();
00412     charge = particle->GetPDGCharge()/eplus;
00413   }
00414   currentKinEnergy = T0;
00415   currentRange = range ;
00416   currentRadLength = couple->GetMaterial()->GetRadlen();
00417 
00418   lambda0 = lambda;
00419   par1 = -1. ;  
00420   par2 = par3 = 0. ;  
00421   tPathLength = truePathLength;
00422 
00423   // this correction needed to run MSC with eIoni and eBrem inactivated
00424   // and makes no harm for a normal run
00425   if(tPathLength > range)
00426     tPathLength = range ;
00427 
00428   G4double tau   = tPathLength/lambda0 ;
00429 
00430   if (tau <= tausmall) return tPathLength;
00431 
00432   G4double zmean = tPathLength;
00433   if (tPathLength < range*dtrl) {
00434     zmean = lambda0*(1.-exp(-tau));
00435     if(tau < taulim) zmean = tPathLength*(1.-0.5*tPathLength/lambda0) ;
00436   } else if(T0 < mass) {
00437     par1 = 1./range ;
00438     par2 = 1./(par1*lambda0) ;
00439     par3 = 1.+par2 ;
00440     zmean = (1.-exp(par3*log(1.-tPathLength/range)))/(par1*par3) ;
00441   } else {
00442     G4LossTableManager* theManager = G4LossTableManager::Instance();
00443     G4double T1 = theManager->GetEnergy(particle,range-tPathLength,couple);
00444     G4double lambda1 ;
00445     if (theLambdaTable) {
00446       G4bool bb;
00447       lambda1 = ((*theLambdaTable)[couple->GetIndex()])->GetValue(T1,bb);
00448     } else {
00449       lambda1 = CrossSection(couple,particle,T1,0.0,1.0);
00450     }
00451     lambda1 = 1.0/lambda1;
00452     par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
00453     par2 = 1./(par1*lambda0) ;
00454     par3 = 1.+par2 ;
00455     zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ;
00456   }
00457 
00458   //  sample z
00459   G4double zPathLength = zmean ;
00460   G4double zt = zmean/tPathLength ;
00461   if (tPathLength >= stepmin && samplez && zt > 0.5 && zt < ztmax)
00462   {
00463     G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
00464     G4double cz1 = 1.+cz ;
00465     G4double u0 = cz/cz1 ;
00466     G4double u,grej ;
00467     do {
00468         u = exp(log(G4UniformRand())/cz1) ;
00469         grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
00470       } while (grej < G4UniformRand()) ;
00471    zPathLength = tPathLength*u ;
00472   }
00473 
00474   return zPathLength ;
00475 }
00476 
00477 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00478 
00479 G4double G4MscModel71::TrueStepLength(G4double geomStepLength)
00480 {
00481   G4double trueLength = geomStepLength;
00482   trueLength = geomStepLength;
00483   if(geomStepLength > lambda0*tausmall)
00484   {
00485     if(par1 <  0.)
00486       trueLength = -lambda0*log(1.-geomStepLength/lambda0) ;
00487     else 
00488     {
00489       if(par1*par3*geomStepLength < 1.)
00490         trueLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
00491       else 
00492         trueLength = currentRange ;
00493     }  
00494   }
00495 
00496   if(trueLength > tPathLength) trueLength = tPathLength;
00497   if(trueLength > currentRange) trueLength = currentRange ;
00498   if(trueLength < geomStepLength) trueLength = geomStepLength;
00499 
00500   return trueLength;
00501 }
00502 
00503 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00504 
00505 void G4MscModel71::SampleSecondaries(std::vector<G4DynamicParticle*>*,
00506                                      const G4MaterialCutsCouple*,
00507                                      const G4DynamicParticle* dynParticle,
00508                                      G4double truestep,
00509                                      G4double safety)
00510 {
00511   G4double kineticEnergy = dynParticle->GetKineticEnergy();
00512   if(kineticEnergy <= 0.0) return;
00513 
00514   G4double cth  = SampleCosineTheta(truestep,kineticEnergy);
00515   G4double sth  = sqrt((1.0 - cth)*(1.0 + cth));
00516   G4double phi  = twopi*G4UniformRand();
00517   G4double dirx = sth*cos(phi);
00518   G4double diry = sth*sin(phi);
00519 
00520   G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
00521   G4ThreeVector newDirection(dirx,diry,cth);
00522   newDirection.rotateUz(oldDirection);
00523   fParticleChange->ProposeMomentumDirection(newDirection);
00524 
00525   /*
00526     const G4ParticleDefinition* pd = dynParticle->GetDefinition();
00527     G4cout << "G4MscModel71: Sample secondary; E(MeV)= " << kineticEnergy/MeV
00528            << " MeV; step(mm)= " << truestep/mm
00529            << ", safety(mm)= " <<  safety/mm << " " << pd->GetParticleName()
00530            << G4endl;
00531   */
00532 
00533   if (latDisplasment && safety > 0.0) {
00534 
00535     G4double r = SampleDisplacement();
00536     if (r > safety) r = safety;
00537 
00538     // sample direction of lateral displacement
00539     G4double phi  = twopi*G4UniformRand();
00540     G4double dirx = std::cos(phi);
00541     G4double diry = std::sin(phi);
00542 
00543     G4ThreeVector newPosition(dirx,diry,0.0);
00544     newPosition.rotateUz(oldDirection);
00545 
00546     // compute new endpoint of the Step
00547     newPosition *= r;
00548     newPosition += *(fParticleChange->GetProposedPosition());
00549 
00550     navigator->LocateGlobalPointWithinVolume(newPosition);
00551 
00552     fParticleChange->ProposePosition(newPosition);
00553   }
00554 }
00555 
00556 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00557 
00558 G4double G4MscModel71::SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
00559 {
00560   G4double cth = 1. ;
00561   G4double tau = trueStepLength/lambda0 ;
00562 
00563   if(trueStepLength >= currentRange*dtrl)
00564     if(par1*trueStepLength < 1.)
00565       tau = -par2*log(1.-par1*trueStepLength) ;
00566     else
00567       tau = taubig ;
00568 
00569   currentTau = tau ;
00570 
00571   if(trueStepLength < stepmin)
00572     cth = exp(-tau) ;
00573   else
00574   {
00575     if (tau >= taubig) cth = -1.+2.*G4UniformRand();
00576     else if (tau >= tausmall)
00577     {
00578         G4double a ;
00579 
00580         // for all particles take the width of the central part
00581         //  from a  parametrization similar to the Highland formula
00582         // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
00583         // here : theta0 = 13.6*MeV*Q*(t/X0)**0.555/(beta*cp) 
00584         const G4double c_highland = 13.6*MeV, corr_highland=0.555 ;
00585         G4double Q = std::abs(charge) ;
00586         G4double xx0 = trueStepLength/currentRadLength;
00587         G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
00588                                KineticEnergy*(KineticEnergy+2.*mass)/
00589                                ((currentKinEnergy+mass)*(KineticEnergy+mass))) ;
00590         G4double theta0 = c_highland*Q*exp(corr_highland*log(xx0))/betacp ;
00591 
00592         if(theta0 > taulim) a = 0.5/(1.-cos(theta0)) ;
00593         else                a = 1.0/(theta0*theta0) ;
00594 
00595         G4double xmeanth = exp(-tau);
00596         G4double xmeanth1 = 1.-xmeanth ;
00597         if(currentTau < taulim) xmeanth1 = tau ;          
00598 
00599         const G4double x1fac1 = exp(-xsi) ;
00600         const G4double x1fac2 = (1.-(1.+xsi)*x1fac1)/(1.-x1fac1) ;
00601         const G4double x1fac3 = 1.3      ; 
00602 
00603         G4double ea,eaa,xmean1 ;
00604         G4double c = 2.,b1 = 2., bx = 2.,
00605                  eb1 = b1, ebx = b1, xmean2 = 0. ;
00606         G4double prob = 1., qprob ;              
00607         G4double x0 = 1.-xsi/a;
00608         G4double oneminusx0=xsi/a ;
00609         G4double oneplusx0=2.+xsi/a ;
00610 
00611         G4double f1x0=1., f2x0=1. ;
00612         const G4double tau0 = 0.10 ;
00613         if(tau > tau0)
00614         {
00615          // 1 model function
00616           a = 1./xmeanth1 ;
00617           ea = exp(-2.*a) ;
00618           eaa= 1.-ea ;
00619           xmean1 = 1.-1./a+2.*ea/eaa ;
00620           prob = 1. ;
00621           qprob = 1. ;
00622         }
00623         else if (x0 <= -1.)
00624         {
00625           // 2 model fuctions only
00626           // in order to have xmean1 > xmeanth -> qprob < 1
00627           x0 = -1.;
00628 
00629           if( a < 1./xmeanth1)
00630             a = 1./xmeanth1 ;
00631           
00632           oneminusx0 = 1.-x0 ;
00633           oneplusx0 =  1.+x0 ;
00634           ea = exp(-a*oneminusx0);
00635           eaa = 1.-ea ;
00636           xmean1 = 1.-1./a+oneminusx0*ea/eaa ;
00637           qprob = xmeanth/xmean1 ;
00638 
00639         }
00640         else
00641         {
00642           // 3 model fuctions
00643           // in order to have xmean1 > xmeanth
00644           if((1.-x1fac2/a) < xmeanth)
00645           {
00646             a = x1fac3*x1fac2/xmeanth1 ;
00647             x0 = 1.-xsi/a ;
00648             oneminusx0=xsi/a ;
00649             oneplusx0=2.-xsi/a ;
00650           }
00651 
00652           ea = x1fac1 ;
00653           eaa = 1.-ea ;
00654           xmean1 = 1.-x1fac2/a ;
00655 
00656           const G4double fctail = factail*1.0 ;
00657           c = 2.+fctail*tau ;
00658           G4double c1 = c-1. ;
00659           G4double c2 = c-2. ;
00660           if(c2 == 0.) c2 = fctail*tausmall ;            
00661 
00662           b = 1.+(c-xsi)/a ;
00663 
00664           b1 = b+1. ;
00665           bx = c/a ;
00666           eb1=exp((c1)*log(b1)) ;
00667           ebx=exp((c1)*log(bx)) ;
00668           xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/c2)/(eb1-ebx) ;
00669 
00670           f1x0 = a*ea/eaa ;
00671           f2x0 = c1*eb1*ebx/(eb1-ebx)/
00672                           exp(c*log(bx)) ;
00673           // from continuity at x=x0
00674           prob = f2x0/(f1x0+f2x0) ;
00675           // from xmean = xmeanth
00676           qprob = (f1x0+f2x0)*xmeanth/(f2x0*xmean1+f1x0*xmean2) ;
00677 
00678         }
00679 
00680         // sampling of costheta
00681         if (G4UniformRand() < qprob)
00682         {
00683           if (G4UniformRand() < prob)
00684              cth = 1.+log(ea+G4UniformRand()*eaa)/a ;
00685           else
00686              cth = b-b1*bx/exp(log(ebx-G4UniformRand()*(ebx-eb1))/(c-1.)) ;
00687         }
00688         else
00689         {
00690           cth = -1.+2.*G4UniformRand();
00691         }
00692     }
00693   }  
00694 
00695   return cth ;
00696 }
00697 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00698 
00699 G4double G4MscModel71::SampleDisplacement()
00700 {
00701   const G4double kappa = 2.5;
00702   const G4double kappapl1 = kappa+1.;
00703   const G4double kappami1 = kappa-1.;
00704   G4double rmean = 0.0;
00705   if (currentTau >= tausmall) {
00706     if (currentTau < taulim) {
00707       rmean = kappa*currentTau*currentTau*currentTau*(1.-kappapl1*currentTau*0.25)/6. ;
00708 
00709     } else {
00710       G4double etau = 0.0;
00711       if (currentTau<taubig) etau = exp(-currentTau);
00712       rmean = -kappa*currentTau;
00713       rmean = -exp(rmean)/(kappa*kappami1);
00714       rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
00715     }
00716     if (rmean>0.) rmean = 2.*lambda0*sqrt(rmean/3.0);
00717     else          rmean = 0.;
00718  }
00719   return rmean;
00720 }
00721 
00722 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00723 
00724 
00725 
00726 

Generated on Fri Apr 11 17:10:10 2008 for Geant4 by  doxygen 1.4.7