| [1197] | 1 | //
|
|---|
| 2 | // ********************************************************************
|
|---|
| 3 | // * License and Disclaimer *
|
|---|
| 4 | // * *
|
|---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of *
|
|---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and *
|
|---|
| 7 | // * conditions of the Geant4 Software License, included in the file *
|
|---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These *
|
|---|
| 9 | // * include a list of copyright holders. *
|
|---|
| 10 | // * *
|
|---|
| 11 | // * Neither the authors of this software system, nor their employing *
|
|---|
| 12 | // * institutes,nor the agencies providing financial support for this *
|
|---|
| 13 | // * work make any representation or warranty, express or implied, *
|
|---|
| 14 | // * regarding this software system or assume any liability for its *
|
|---|
| 15 | // * use. Please see the license in the file LICENSE and URL above *
|
|---|
| 16 | // * for the full disclaimer and the limitation of liability. *
|
|---|
| 17 | // * *
|
|---|
| 18 | // * This code implementation is the result of the scientific and *
|
|---|
| 19 | // * technical work of the GEANT4 collaboration. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| 26 | //
|
|---|
| 27 | // $Id: G4UrbanMscModel93.cc,v 1.1 2009/11/01 13:05:01 vnivanch Exp $
|
|---|
| 28 | // GEANT4 tag $Name: geant4-09-03-cand-01 $
|
|---|
| 29 | //
|
|---|
| 30 | // -------------------------------------------------------------------
|
|---|
| 31 | //
|
|---|
| 32 | // GEANT4 Class file
|
|---|
| 33 | //
|
|---|
| 34 | //
|
|---|
| 35 | // File name: G4UrbanMscModel93
|
|---|
| 36 | //
|
|---|
| 37 | // Author: Laszlo Urban
|
|---|
| 38 | //
|
|---|
| 39 | // Creation date: 06.03.2008
|
|---|
| 40 | //
|
|---|
| 41 | // Modifications:
|
|---|
| 42 | //
|
|---|
| 43 | // 06-03-2008 starting point : G4UrbanMscModel2 = G4UrbanMscModel 9.1 ref 02
|
|---|
| 44 | //
|
|---|
| 45 | // 13-03-08 Bug in SampleScattering (which caused lateral asymmetry) fixed
|
|---|
| 46 | // (L.Urban)
|
|---|
| 47 | //
|
|---|
| 48 | // 14-03-08 Simplification of step limitation in ComputeTruePathLengthLimit,
|
|---|
| 49 | // + tlimitmin is the same for UseDistancetoBoundary and
|
|---|
| 50 | // UseSafety (L.Urban)
|
|---|
| 51 | //
|
|---|
| 52 | // 16-03-08 Reorganization of SampleCosineTheta + new method SimpleScattering
|
|---|
| 53 | // SimpleScattering is used if the relative energy loss is too big
|
|---|
| 54 | // or theta0 is too big (see data members rellossmax, theta0max)
|
|---|
| 55 | // (L.Urban)
|
|---|
| 56 | //
|
|---|
| 57 | // 17-03-08 tuning of the correction factor in ComputeTheta0 (L.Urban)
|
|---|
| 58 | //
|
|---|
| 59 | // 19-03-08 exponent c of the 'tail' model function is not equal to 2 any more,
|
|---|
| 60 | // value of c has been extracted from some e- scattering data (L.Urban)
|
|---|
| 61 | //
|
|---|
| 62 | // 24-03-08 Step limitation in ComputeTruePathLengthLimit has been
|
|---|
| 63 | // simplified further + some data members have been removed (L.Urban)
|
|---|
| 64 | //
|
|---|
| 65 | // 24-07-08 central part of scattering angle (theta0) has been tuned
|
|---|
| 66 | // tail of the scattering angle distribution has been tuned
|
|---|
| 67 | // using some e- and proton scattering data
|
|---|
| 68 | //
|
|---|
| 69 | // 05-08-08 bugfix in ComputeTruePathLengthLimit (L.Urban)
|
|---|
| 70 | //
|
|---|
| 71 | // 09-10-08 theta0 and tail have been retuned using some e-,mu,proton
|
|---|
| 72 | // scattering data (L.Urban)
|
|---|
| 73 | // + single scattering without path length correction for
|
|---|
| 74 | // small steps (t < tlimitmin, for UseDistanceToBoundary only)
|
|---|
| 75 | //
|
|---|
| 76 | // 15-10-08 Moliere-Bethe screening in the single scattering part(L.Urban)
|
|---|
| 77 | //
|
|---|
| 78 | // 17-10-08 stepping similar to that in model (9.1) for UseSafety case
|
|---|
| 79 | // for e+/e- in order to speed up the code for calorimeters
|
|---|
| 80 | //
|
|---|
| 81 | // 23-10-08 bugfix in the screeningparameter of the single scattering part,
|
|---|
| 82 | // some technical change in order to speed up the code (UpdateCache)
|
|---|
| 83 | //
|
|---|
| 84 | // 27-10-08 bugfix in ComputeTruePathLengthLimit (affects UseDistanceToBoundary
|
|---|
| 85 | // stepping type only) (L.Urban)
|
|---|
| 86 | //
|
|---|
| 87 | // 28-10-09 V.Ivanchenko moved G4UrbanMscModel2 to G4UrbanMscModel93,
|
|---|
| 88 | // now it is a frozen version of the Urban model corresponding
|
|---|
| 89 | // to g4 9.3
|
|---|
| 90 |
|
|---|
| 91 | // Class Description:
|
|---|
| 92 | //
|
|---|
| 93 | // Implementation of the model of multiple scattering based on
|
|---|
| 94 | // H.W.Lewis Phys Rev 78 (1950) 526 and others
|
|---|
| 95 |
|
|---|
| 96 | // -------------------------------------------------------------------
|
|---|
| 97 | // In its present form the model can be used for simulation
|
|---|
| 98 | // of the e-/e+, muon and charged hadron multiple scattering
|
|---|
| 99 | //
|
|---|
| 100 |
|
|---|
| 101 |
|
|---|
| 102 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 103 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 104 |
|
|---|
| 105 | #include "G4UrbanMscModel93.hh"
|
|---|
| 106 | #include "Randomize.hh"
|
|---|
| 107 | #include "G4Electron.hh"
|
|---|
| 108 | #include "G4LossTableManager.hh"
|
|---|
| 109 | #include "G4ParticleChangeForMSC.hh"
|
|---|
| 110 |
|
|---|
| 111 | #include "G4Poisson.hh"
|
|---|
| 112 | #include "globals.hh"
|
|---|
| 113 |
|
|---|
| 114 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 115 |
|
|---|
| 116 | using namespace std;
|
|---|
| 117 |
|
|---|
| 118 | G4UrbanMscModel93::G4UrbanMscModel93(const G4String& nam)
|
|---|
| 119 | : G4VMscModel(nam),
|
|---|
| 120 | isInitialized(false)
|
|---|
| 121 | {
|
|---|
| 122 | masslimite = 0.6*MeV;
|
|---|
| 123 | lambdalimit = 1.*mm;
|
|---|
| 124 | fr = 0.02;
|
|---|
| 125 | //facsafety = 0.3;
|
|---|
| 126 | taubig = 8.0;
|
|---|
| 127 | tausmall = 1.e-16;
|
|---|
| 128 | taulim = 1.e-6;
|
|---|
| 129 | currentTau = taulim;
|
|---|
| 130 | tlimitminfix = 1.e-6*mm;
|
|---|
| 131 | stepmin = tlimitminfix;
|
|---|
| 132 | smallstep = 1.e10;
|
|---|
| 133 | currentRange = 0. ;
|
|---|
| 134 | rangeinit = 0.;
|
|---|
| 135 | tlimit = 1.e10*mm;
|
|---|
| 136 | tlimitmin = 10.*tlimitminfix;
|
|---|
| 137 | tgeom = 1.e50*mm;
|
|---|
| 138 | geombig = 1.e50*mm;
|
|---|
| 139 | geommin = 1.e-3*mm;
|
|---|
| 140 | geomlimit = geombig;
|
|---|
| 141 | presafety = 0.*mm;
|
|---|
| 142 |
|
|---|
| 143 | y = 0.;
|
|---|
| 144 |
|
|---|
| 145 | Zold = 0.;
|
|---|
| 146 | Zeff = 1.;
|
|---|
| 147 | Z2 = 1.;
|
|---|
| 148 | Z23 = 1.;
|
|---|
| 149 | lnZ = 0.;
|
|---|
| 150 | coeffth1 = 0.;
|
|---|
| 151 | coeffth2 = 0.;
|
|---|
| 152 | coeffc1 = 0.;
|
|---|
| 153 | coeffc2 = 0.;
|
|---|
| 154 | scr1ini = fine_structure_const*fine_structure_const*
|
|---|
| 155 | electron_mass_c2*electron_mass_c2/(0.885*0.885*4.*pi);
|
|---|
| 156 | scr2ini = 3.76*fine_structure_const*fine_structure_const;
|
|---|
| 157 | scr1 = 0.;
|
|---|
| 158 | scr2 = 0.;
|
|---|
| 159 |
|
|---|
| 160 | theta0max = pi/6.;
|
|---|
| 161 | rellossmax = 0.50;
|
|---|
| 162 | third = 1./3.;
|
|---|
| 163 | particle = 0;
|
|---|
| 164 | theManager = G4LossTableManager::Instance();
|
|---|
| 165 | inside = false;
|
|---|
| 166 | insideskin = false;
|
|---|
| 167 |
|
|---|
| 168 | }
|
|---|
| 169 |
|
|---|
| 170 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 171 |
|
|---|
| 172 | G4UrbanMscModel93::~G4UrbanMscModel93()
|
|---|
| 173 | {}
|
|---|
| 174 |
|
|---|
| 175 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 176 |
|
|---|
| 177 | void G4UrbanMscModel93::Initialise(const G4ParticleDefinition* p,
|
|---|
| 178 | const G4DataVector&)
|
|---|
| 179 | {
|
|---|
| 180 | skindepth = skin*stepmin;
|
|---|
| 181 | if(isInitialized) return;
|
|---|
| 182 | // set values of some data members
|
|---|
| 183 | SetParticle(p);
|
|---|
| 184 |
|
|---|
| 185 | fParticleChange = GetParticleChangeForMSC();
|
|---|
| 186 | InitialiseSafetyHelper();
|
|---|
| 187 |
|
|---|
| 188 | isInitialized = true;
|
|---|
| 189 | }
|
|---|
| 190 |
|
|---|
| 191 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 192 |
|
|---|
| 193 | G4double G4UrbanMscModel93::ComputeCrossSectionPerAtom(
|
|---|
| 194 | const G4ParticleDefinition* part,
|
|---|
| 195 | G4double KineticEnergy,
|
|---|
| 196 | G4double AtomicNumber,G4double,
|
|---|
| 197 | G4double, G4double)
|
|---|
| 198 | {
|
|---|
| 199 | const G4double sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
|
|---|
| 200 | const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
|
|---|
| 201 | Bohr_radius*Bohr_radius/(hbarc*hbarc);
|
|---|
| 202 | const G4double epsmin = 1.e-4 , epsmax = 1.e10;
|
|---|
| 203 |
|
|---|
| 204 | const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38., 47.,
|
|---|
| 205 | 50., 56., 64., 74., 79., 82. };
|
|---|
| 206 |
|
|---|
| 207 | const G4double Tdat[22] = { 100*eV, 200*eV, 400*eV, 700*eV,
|
|---|
| 208 | 1*keV, 2*keV, 4*keV, 7*keV,
|
|---|
| 209 | 10*keV, 20*keV, 40*keV, 70*keV,
|
|---|
| 210 | 100*keV, 200*keV, 400*keV, 700*keV,
|
|---|
| 211 | 1*MeV, 2*MeV, 4*MeV, 7*MeV,
|
|---|
| 212 | 10*MeV, 20*MeV};
|
|---|
| 213 |
|
|---|
| 214 | // corr. factors for e-/e+ lambda for T <= Tlim
|
|---|
| 215 | G4double celectron[15][22] =
|
|---|
| 216 | {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
|
|---|
| 217 | 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
|
|---|
| 218 | 1.112,1.108,1.100,1.093,1.089,1.087 },
|
|---|
| 219 | {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
|
|---|
| 220 | 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
|
|---|
| 221 | 1.109,1.105,1.097,1.090,1.086,1.082 },
|
|---|
| 222 | {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
|
|---|
| 223 | 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
|
|---|
| 224 | 1.131,1.124,1.113,1.104,1.099,1.098 },
|
|---|
| 225 | {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
|
|---|
| 226 | 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
|
|---|
| 227 | 1.112,1.105,1.096,1.089,1.085,1.098 },
|
|---|
| 228 | {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
|
|---|
| 229 | 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
|
|---|
| 230 | 1.073,1.070,1.064,1.059,1.056,1.056 },
|
|---|
| 231 | {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
|
|---|
| 232 | 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
|
|---|
| 233 | 1.074,1.070,1.063,1.059,1.056,1.052 },
|
|---|
| 234 | {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
|
|---|
| 235 | 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
|
|---|
| 236 | 1.068,1.064,1.059,1.054,1.051,1.050 },
|
|---|
| 237 | {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
|
|---|
| 238 | 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
|
|---|
| 239 | 1.039,1.037,1.034,1.031,1.030,1.036 },
|
|---|
| 240 | {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
|
|---|
| 241 | 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
|
|---|
| 242 | 1.031,1.028,1.024,1.022,1.021,1.024 },
|
|---|
| 243 | {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
|
|---|
| 244 | 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
|
|---|
| 245 | 1.020,1.017,1.015,1.013,1.013,1.020 },
|
|---|
| 246 | {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
|
|---|
| 247 | 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
|
|---|
| 248 | 0.995,0.993,0.993,0.993,0.993,1.011 },
|
|---|
| 249 | {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
|
|---|
| 250 | 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
|
|---|
| 251 | 0.974,0.972,0.973,0.974,0.975,0.987 },
|
|---|
| 252 | {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
|
|---|
| 253 | 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
|
|---|
| 254 | 0.950,0.947,0.949,0.952,0.954,0.963 },
|
|---|
| 255 | {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
|
|---|
| 256 | 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
|
|---|
| 257 | 0.941,0.938,0.940,0.944,0.946,0.954 },
|
|---|
| 258 | {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
|
|---|
| 259 | 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
|
|---|
| 260 | 0.933,0.930,0.933,0.936,0.939,0.949 }};
|
|---|
| 261 |
|
|---|
| 262 | G4double cpositron[15][22] = {
|
|---|
| 263 | {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
|
|---|
| 264 | 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
|
|---|
| 265 | 1.131,1.126,1.117,1.108,1.103,1.100 },
|
|---|
| 266 | {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
|
|---|
| 267 | 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
|
|---|
| 268 | 1.138,1.132,1.122,1.113,1.108,1.102 },
|
|---|
| 269 | {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
|
|---|
| 270 | 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
|
|---|
| 271 | 1.203,1.190,1.173,1.159,1.151,1.145 },
|
|---|
| 272 | {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
|
|---|
| 273 | 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
|
|---|
| 274 | 1.225,1.210,1.191,1.175,1.166,1.174 },
|
|---|
| 275 | {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
|
|---|
| 276 | 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
|
|---|
| 277 | 1.217,1.203,1.184,1.169,1.160,1.151 },
|
|---|
| 278 | {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
|
|---|
| 279 | 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
|
|---|
| 280 | 1.237,1.222,1.201,1.184,1.174,1.159 },
|
|---|
| 281 | {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
|
|---|
| 282 | 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
|
|---|
| 283 | 1.252,1.234,1.212,1.194,1.183,1.170 },
|
|---|
| 284 | {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
|
|---|
| 285 | 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
|
|---|
| 286 | 1.254,1.237,1.214,1.195,1.185,1.179 },
|
|---|
| 287 | {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
|
|---|
| 288 | 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
|
|---|
| 289 | 1.312,1.288,1.258,1.235,1.221,1.205 },
|
|---|
| 290 | {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
|
|---|
| 291 | 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
|
|---|
| 292 | 1.320,1.294,1.264,1.240,1.226,1.214 },
|
|---|
| 293 | {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
|
|---|
| 294 | 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
|
|---|
| 295 | 1.328,1.302,1.270,1.245,1.231,1.233 },
|
|---|
| 296 | {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
|
|---|
| 297 | 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
|
|---|
| 298 | 1.361,1.330,1.294,1.267,1.251,1.239 },
|
|---|
| 299 | {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
|
|---|
| 300 | 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
|
|---|
| 301 | 1.409,1.372,1.330,1.298,1.280,1.258 },
|
|---|
| 302 | {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
|
|---|
| 303 | 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
|
|---|
| 304 | 1.442,1.400,1.354,1.319,1.299,1.272 },
|
|---|
| 305 | {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
|
|---|
| 306 | 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
|
|---|
| 307 | 1.456,1.412,1.364,1.328,1.307,1.282 }};
|
|---|
| 308 |
|
|---|
| 309 | //data/corrections for T > Tlim
|
|---|
| 310 | G4double Tlim = 10.*MeV;
|
|---|
| 311 | G4double beta2lim = Tlim*(Tlim+2.*electron_mass_c2)/
|
|---|
| 312 | ((Tlim+electron_mass_c2)*(Tlim+electron_mass_c2));
|
|---|
| 313 | G4double bg2lim = Tlim*(Tlim+2.*electron_mass_c2)/
|
|---|
| 314 | (electron_mass_c2*electron_mass_c2);
|
|---|
| 315 |
|
|---|
| 316 | G4double sig0[15] = {0.2672*barn, 0.5922*barn, 2.653*barn, 6.235*barn,
|
|---|
| 317 | 11.69*barn , 13.24*barn , 16.12*barn, 23.00*barn ,
|
|---|
| 318 | 35.13*barn , 39.95*barn , 50.85*barn, 67.19*barn ,
|
|---|
| 319 | 91.15*barn , 104.4*barn , 113.1*barn};
|
|---|
| 320 |
|
|---|
| 321 | G4double hecorr[15] = {120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
|
|---|
| 322 | 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
|
|---|
| 323 | -22.30};
|
|---|
| 324 |
|
|---|
| 325 | G4double sigma;
|
|---|
| 326 | SetParticle(part);
|
|---|
| 327 |
|
|---|
| 328 | G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
|
|---|
| 329 |
|
|---|
| 330 | // correction if particle .ne. e-/e+
|
|---|
| 331 | // compute equivalent kinetic energy
|
|---|
| 332 | // lambda depends on p*beta ....
|
|---|
| 333 |
|
|---|
| 334 | G4double eKineticEnergy = KineticEnergy;
|
|---|
| 335 |
|
|---|
| 336 | if(mass > electron_mass_c2)
|
|---|
| 337 | {
|
|---|
| 338 | G4double TAU = KineticEnergy/mass ;
|
|---|
| 339 | G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
|
|---|
| 340 | G4double w = c-2. ;
|
|---|
| 341 | G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
|
|---|
| 342 | eKineticEnergy = electron_mass_c2*tau ;
|
|---|
| 343 | }
|
|---|
| 344 |
|
|---|
| 345 | G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
|
|---|
| 346 | G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
|
|---|
| 347 | /(eTotalEnergy*eTotalEnergy);
|
|---|
| 348 | G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
|
|---|
| 349 | /(electron_mass_c2*electron_mass_c2);
|
|---|
| 350 |
|
|---|
| 351 | G4double eps = epsfactor*bg2/Z23;
|
|---|
| 352 |
|
|---|
| 353 | if (eps<epsmin) sigma = 2.*eps*eps;
|
|---|
| 354 | else if(eps<epsmax) sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
|
|---|
| 355 | else sigma = log(2.*eps)-1.+1./eps;
|
|---|
| 356 |
|
|---|
| 357 | sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
|
|---|
| 358 |
|
|---|
| 359 | // interpolate in AtomicNumber and beta2
|
|---|
| 360 | G4double c1,c2,cc1,cc2,corr;
|
|---|
| 361 |
|
|---|
| 362 | // get bin number in Z
|
|---|
| 363 | G4int iZ = 14;
|
|---|
| 364 | while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
|
|---|
| 365 | if (iZ==14) iZ = 13;
|
|---|
| 366 | if (iZ==-1) iZ = 0 ;
|
|---|
| 367 |
|
|---|
| 368 | G4double Z1 = Zdat[iZ];
|
|---|
| 369 | G4double Z2 = Zdat[iZ+1];
|
|---|
| 370 | G4double ratZ = (AtomicNumber-Z1)*(AtomicNumber+Z1)/
|
|---|
| 371 | ((Z2-Z1)*(Z2+Z1));
|
|---|
| 372 |
|
|---|
| 373 | if(eKineticEnergy <= Tlim)
|
|---|
| 374 | {
|
|---|
| 375 | // get bin number in T (beta2)
|
|---|
| 376 | G4int iT = 21;
|
|---|
| 377 | while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
|
|---|
| 378 | if(iT==21) iT = 20;
|
|---|
| 379 | if(iT==-1) iT = 0 ;
|
|---|
| 380 |
|
|---|
| 381 | // calculate betasquare values
|
|---|
| 382 | G4double T = Tdat[iT], E = T + electron_mass_c2;
|
|---|
| 383 | G4double b2small = T*(E+electron_mass_c2)/(E*E);
|
|---|
| 384 |
|
|---|
| 385 | T = Tdat[iT+1]; E = T + electron_mass_c2;
|
|---|
| 386 | G4double b2big = T*(E+electron_mass_c2)/(E*E);
|
|---|
| 387 | G4double ratb2 = (beta2-b2small)/(b2big-b2small);
|
|---|
| 388 |
|
|---|
| 389 | if (charge < 0.)
|
|---|
| 390 | {
|
|---|
| 391 | c1 = celectron[iZ][iT];
|
|---|
| 392 | c2 = celectron[iZ+1][iT];
|
|---|
| 393 | cc1 = c1+ratZ*(c2-c1);
|
|---|
| 394 |
|
|---|
| 395 | c1 = celectron[iZ][iT+1];
|
|---|
| 396 | c2 = celectron[iZ+1][iT+1];
|
|---|
| 397 | cc2 = c1+ratZ*(c2-c1);
|
|---|
| 398 |
|
|---|
| 399 | corr = cc1+ratb2*(cc2-cc1);
|
|---|
| 400 |
|
|---|
| 401 | sigma *= sigmafactor/corr;
|
|---|
| 402 | }
|
|---|
| 403 | else
|
|---|
| 404 | {
|
|---|
| 405 | c1 = cpositron[iZ][iT];
|
|---|
| 406 | c2 = cpositron[iZ+1][iT];
|
|---|
| 407 | cc1 = c1+ratZ*(c2-c1);
|
|---|
| 408 |
|
|---|
| 409 | c1 = cpositron[iZ][iT+1];
|
|---|
| 410 | c2 = cpositron[iZ+1][iT+1];
|
|---|
| 411 | cc2 = c1+ratZ*(c2-c1);
|
|---|
| 412 |
|
|---|
| 413 | corr = cc1+ratb2*(cc2-cc1);
|
|---|
| 414 |
|
|---|
| 415 | sigma *= sigmafactor/corr;
|
|---|
| 416 | }
|
|---|
| 417 | }
|
|---|
| 418 | else
|
|---|
| 419 | {
|
|---|
| 420 | c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
|
|---|
| 421 | c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
|
|---|
| 422 | if((AtomicNumber >= Z1) && (AtomicNumber <= Z2))
|
|---|
| 423 | sigma = c1+ratZ*(c2-c1) ;
|
|---|
| 424 | else if(AtomicNumber < Z1)
|
|---|
| 425 | sigma = AtomicNumber*AtomicNumber*c1/(Z1*Z1);
|
|---|
| 426 | else if(AtomicNumber > Z2)
|
|---|
| 427 | sigma = AtomicNumber*AtomicNumber*c2/(Z2*Z2);
|
|---|
| 428 | }
|
|---|
| 429 | return sigma;
|
|---|
| 430 |
|
|---|
| 431 | }
|
|---|
| 432 |
|
|---|
| 433 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 434 |
|
|---|
| 435 | G4double G4UrbanMscModel93::ComputeTruePathLengthLimit(
|
|---|
| 436 | const G4Track& track,
|
|---|
| 437 | G4PhysicsTable* theTable,
|
|---|
| 438 | G4double currentMinimalStep)
|
|---|
| 439 | {
|
|---|
| 440 | tPathLength = currentMinimalStep;
|
|---|
| 441 | G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
|
|---|
| 442 | G4StepStatus stepStatus = sp->GetStepStatus();
|
|---|
| 443 |
|
|---|
| 444 | const G4DynamicParticle* dp = track.GetDynamicParticle();
|
|---|
| 445 |
|
|---|
| 446 | if(stepStatus == fUndefined) {
|
|---|
| 447 | inside = false;
|
|---|
| 448 | insideskin = false;
|
|---|
| 449 | tlimit = geombig;
|
|---|
| 450 | SetParticle( dp->GetDefinition() );
|
|---|
| 451 | }
|
|---|
| 452 |
|
|---|
| 453 | theLambdaTable = theTable;
|
|---|
| 454 | couple = track.GetMaterialCutsCouple();
|
|---|
| 455 | currentMaterialIndex = couple->GetIndex();
|
|---|
| 456 | currentKinEnergy = dp->GetKineticEnergy();
|
|---|
| 457 | currentRange =
|
|---|
| 458 | theManager->GetRangeFromRestricteDEDX(particle,currentKinEnergy,couple);
|
|---|
| 459 | lambda0 = GetLambda(currentKinEnergy);
|
|---|
| 460 |
|
|---|
| 461 | // stop here if small range particle
|
|---|
| 462 | if(inside) return tPathLength;
|
|---|
| 463 |
|
|---|
| 464 | if(tPathLength > currentRange) tPathLength = currentRange;
|
|---|
| 465 |
|
|---|
| 466 | presafety = sp->GetSafety();
|
|---|
| 467 |
|
|---|
| 468 | //G4cout << "G4Urban2::StepLimit tPathLength= "
|
|---|
| 469 | // <<tPathLength<<" safety= " << presafety
|
|---|
| 470 | // << " range= " <<currentRange<< " lambda= "<<lambda0
|
|---|
| 471 | // << " Alg: " << steppingAlgorithm <<G4endl;
|
|---|
| 472 |
|
|---|
| 473 | // far from geometry boundary
|
|---|
| 474 | if(currentRange < presafety)
|
|---|
| 475 | {
|
|---|
| 476 | inside = true;
|
|---|
| 477 | return tPathLength;
|
|---|
| 478 | }
|
|---|
| 479 |
|
|---|
| 480 | // standard version
|
|---|
| 481 | //
|
|---|
| 482 | if (steppingAlgorithm == fUseDistanceToBoundary)
|
|---|
| 483 | {
|
|---|
| 484 | //compute geomlimit and presafety
|
|---|
| 485 | G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
|
|---|
| 486 |
|
|---|
| 487 | // is it far from boundary ?
|
|---|
| 488 | if(currentRange < presafety)
|
|---|
| 489 | {
|
|---|
| 490 | inside = true;
|
|---|
| 491 | return tPathLength;
|
|---|
| 492 | }
|
|---|
| 493 |
|
|---|
| 494 | smallstep += 1.;
|
|---|
| 495 | insideskin = false;
|
|---|
| 496 |
|
|---|
| 497 | if((stepStatus == fGeomBoundary) || (stepStatus == fUndefined))
|
|---|
| 498 | {
|
|---|
| 499 | rangeinit = currentRange;
|
|---|
| 500 | if(stepStatus == fUndefined) smallstep = 1.e10;
|
|---|
| 501 | else smallstep = 1.;
|
|---|
| 502 |
|
|---|
| 503 | //define stepmin here (it depends on lambda!)
|
|---|
| 504 | //rough estimation of lambda_elastic/lambda_transport
|
|---|
| 505 | G4double rat = currentKinEnergy/MeV ;
|
|---|
| 506 | rat = 1.e-3/(rat*(10.+rat)) ;
|
|---|
| 507 | //stepmin ~ lambda_elastic
|
|---|
| 508 | stepmin = rat*lambda0;
|
|---|
| 509 | skindepth = skin*stepmin;
|
|---|
| 510 | //define tlimitmin
|
|---|
| 511 | tlimitmin = 10.*stepmin;
|
|---|
| 512 | if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
|
|---|
| 513 | //G4cout << "rangeinit= " << rangeinit << " stepmin= " << stepmin
|
|---|
| 514 | // << " tlimitmin= " << tlimitmin << " geomlimit= " << geomlimit <<G4endl;
|
|---|
| 515 | // constraint from the geometry
|
|---|
| 516 | if((geomlimit < geombig) && (geomlimit > geommin))
|
|---|
| 517 | {
|
|---|
| 518 | // geomlimit is a geometrical step length
|
|---|
| 519 | // transform it to true path length (estimation)
|
|---|
| 520 | if((1.-geomlimit/lambda0) > 0.)
|
|---|
| 521 | geomlimit = -lambda0*log(1.-geomlimit/lambda0)+tlimitmin ;
|
|---|
| 522 |
|
|---|
| 523 | if(stepStatus == fGeomBoundary)
|
|---|
| 524 | tgeom = geomlimit/facgeom;
|
|---|
| 525 | else
|
|---|
| 526 | tgeom = 2.*geomlimit/facgeom;
|
|---|
| 527 | }
|
|---|
| 528 | else
|
|---|
| 529 | tgeom = geombig;
|
|---|
| 530 | }
|
|---|
| 531 |
|
|---|
| 532 |
|
|---|
| 533 | //step limit
|
|---|
| 534 | tlimit = facrange*rangeinit;
|
|---|
| 535 | if(tlimit < facsafety*presafety)
|
|---|
| 536 | tlimit = facsafety*presafety;
|
|---|
| 537 |
|
|---|
| 538 | //lower limit for tlimit
|
|---|
| 539 | if(tlimit < tlimitmin) tlimit = tlimitmin;
|
|---|
| 540 |
|
|---|
| 541 | if(tlimit > tgeom) tlimit = tgeom;
|
|---|
| 542 |
|
|---|
| 543 | //G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit
|
|---|
| 544 | // << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
|
|---|
| 545 |
|
|---|
| 546 | // shortcut
|
|---|
| 547 | if((tPathLength < tlimit) && (tPathLength < presafety) &&
|
|---|
| 548 | (smallstep >= skin) && (tPathLength < geomlimit-0.999*skindepth))
|
|---|
| 549 | return tPathLength;
|
|---|
| 550 |
|
|---|
| 551 | // step reduction near to boundary
|
|---|
| 552 | if(smallstep < skin)
|
|---|
| 553 | {
|
|---|
| 554 | tlimit = stepmin;
|
|---|
| 555 | insideskin = true;
|
|---|
| 556 | }
|
|---|
| 557 | else if(geomlimit < geombig)
|
|---|
| 558 | {
|
|---|
| 559 | if(geomlimit > skindepth)
|
|---|
| 560 | {
|
|---|
| 561 | if(tlimit > geomlimit-0.999*skindepth)
|
|---|
| 562 | tlimit = geomlimit-0.999*skindepth;
|
|---|
| 563 | }
|
|---|
| 564 | else
|
|---|
| 565 | {
|
|---|
| 566 | insideskin = true;
|
|---|
| 567 | if(tlimit > stepmin) tlimit = stepmin;
|
|---|
| 568 | }
|
|---|
| 569 | }
|
|---|
| 570 |
|
|---|
| 571 | if(tlimit < stepmin) tlimit = stepmin;
|
|---|
| 572 |
|
|---|
| 573 | if(tPathLength > tlimit) tPathLength = tlimit ;
|
|---|
| 574 |
|
|---|
| 575 | }
|
|---|
| 576 | // for 'normal' simulation with or without magnetic field
|
|---|
| 577 | // there no small step/single scattering at boundaries
|
|---|
| 578 | else if(steppingAlgorithm == fUseSafety)
|
|---|
| 579 | {
|
|---|
| 580 | // compute presafety again if presafety <= 0 and no boundary
|
|---|
| 581 | // i.e. when it is needed for optimization purposes
|
|---|
| 582 | if((stepStatus != fGeomBoundary) && (presafety < tlimitminfix))
|
|---|
| 583 | presafety = ComputeSafety(sp->GetPosition(),tPathLength);
|
|---|
| 584 |
|
|---|
| 585 | // is far from boundary
|
|---|
| 586 | if(currentRange < presafety)
|
|---|
| 587 | {
|
|---|
| 588 | inside = true;
|
|---|
| 589 | return tPathLength;
|
|---|
| 590 | }
|
|---|
| 591 |
|
|---|
| 592 | if((stepStatus == fGeomBoundary) || (stepStatus == fUndefined))
|
|---|
| 593 | {
|
|---|
| 594 | rangeinit = currentRange;
|
|---|
| 595 | fr = facrange;
|
|---|
| 596 | // 9.1 like stepping for e+/e- only (not for muons,hadrons)
|
|---|
| 597 | if(mass < masslimite)
|
|---|
| 598 | {
|
|---|
| 599 | if(lambda0 > currentRange)
|
|---|
| 600 | rangeinit = lambda0;
|
|---|
| 601 | if(lambda0 > lambdalimit)
|
|---|
| 602 | fr *= 0.75+0.25*lambda0/lambdalimit;
|
|---|
| 603 | }
|
|---|
| 604 |
|
|---|
| 605 | //lower limit for tlimit
|
|---|
| 606 | G4double rat = currentKinEnergy/MeV ;
|
|---|
| 607 | rat = 1.e-3/(rat*(10.+rat)) ;
|
|---|
| 608 | tlimitmin = 10.*lambda0*rat;
|
|---|
| 609 | if(tlimitmin < tlimitminfix) tlimitmin = tlimitminfix;
|
|---|
| 610 | }
|
|---|
| 611 | //step limit
|
|---|
| 612 | tlimit = fr*rangeinit;
|
|---|
| 613 |
|
|---|
| 614 | if(tlimit < facsafety*presafety)
|
|---|
| 615 | tlimit = facsafety*presafety;
|
|---|
| 616 |
|
|---|
| 617 | //lower limit for tlimit
|
|---|
| 618 | if(tlimit < tlimitmin) tlimit = tlimitmin;
|
|---|
| 619 |
|
|---|
| 620 | if(tPathLength > tlimit) tPathLength = tlimit;
|
|---|
| 621 | }
|
|---|
| 622 |
|
|---|
| 623 | // version similar to 7.1 (needed for some experiments)
|
|---|
| 624 | else
|
|---|
| 625 | {
|
|---|
| 626 | if (stepStatus == fGeomBoundary)
|
|---|
| 627 | {
|
|---|
| 628 | if (currentRange > lambda0) tlimit = facrange*currentRange;
|
|---|
| 629 | else tlimit = facrange*lambda0;
|
|---|
| 630 |
|
|---|
| 631 | if(tlimit < tlimitmin) tlimit = tlimitmin;
|
|---|
| 632 | if(tPathLength > tlimit) tPathLength = tlimit;
|
|---|
| 633 | }
|
|---|
| 634 | }
|
|---|
| 635 | //G4cout << "tPathLength= " << tPathLength
|
|---|
| 636 | // << " currentMinimalStep= " << currentMinimalStep << G4endl;
|
|---|
| 637 | return tPathLength ;
|
|---|
| 638 | }
|
|---|
| 639 |
|
|---|
| 640 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 641 |
|
|---|
| 642 | G4double G4UrbanMscModel93::ComputeGeomPathLength(G4double)
|
|---|
| 643 | {
|
|---|
| 644 | lambdaeff = lambda0;
|
|---|
| 645 | par1 = -1. ;
|
|---|
| 646 | par2 = par3 = 0. ;
|
|---|
| 647 |
|
|---|
| 648 | // do the true -> geom transformation
|
|---|
| 649 | zPathLength = tPathLength;
|
|---|
| 650 |
|
|---|
| 651 | // z = t for very small tPathLength
|
|---|
| 652 | if(tPathLength < tlimitminfix) return zPathLength;
|
|---|
| 653 |
|
|---|
| 654 | // this correction needed to run MSC with eIoni and eBrem inactivated
|
|---|
| 655 | // and makes no harm for a normal run
|
|---|
| 656 | if(tPathLength > currentRange)
|
|---|
| 657 | tPathLength = currentRange ;
|
|---|
| 658 |
|
|---|
| 659 | G4double tau = tPathLength/lambda0 ;
|
|---|
| 660 |
|
|---|
| 661 | if ((tau <= tausmall) || insideskin) {
|
|---|
| 662 | zPathLength = tPathLength;
|
|---|
| 663 | if(zPathLength > lambda0) zPathLength = lambda0;
|
|---|
| 664 | return zPathLength;
|
|---|
| 665 | }
|
|---|
| 666 |
|
|---|
| 667 | G4double zmean = tPathLength;
|
|---|
| 668 | if (tPathLength < currentRange*dtrl) {
|
|---|
| 669 | if(tau < taulim) zmean = tPathLength*(1.-0.5*tau) ;
|
|---|
| 670 | else zmean = lambda0*(1.-exp(-tau));
|
|---|
| 671 | } else if(currentKinEnergy < mass) {
|
|---|
| 672 | par1 = 1./currentRange ;
|
|---|
| 673 | par2 = 1./(par1*lambda0) ;
|
|---|
| 674 | par3 = 1.+par2 ;
|
|---|
| 675 | if(tPathLength < currentRange)
|
|---|
| 676 | zmean = (1.-exp(par3*log(1.-tPathLength/currentRange)))/(par1*par3) ;
|
|---|
| 677 | else
|
|---|
| 678 | zmean = 1./(par1*par3) ;
|
|---|
| 679 | } else {
|
|---|
| 680 | G4double T1 = theManager->GetEnergy(particle,currentRange-tPathLength,couple);
|
|---|
| 681 | G4double lambda1 = GetLambda(T1);
|
|---|
| 682 |
|
|---|
| 683 | par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
|
|---|
| 684 | par2 = 1./(par1*lambda0) ;
|
|---|
| 685 | par3 = 1.+par2 ;
|
|---|
| 686 | zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ;
|
|---|
| 687 | }
|
|---|
| 688 |
|
|---|
| 689 | zPathLength = zmean ;
|
|---|
| 690 |
|
|---|
| 691 | // sample z
|
|---|
| 692 | if(samplez)
|
|---|
| 693 | {
|
|---|
| 694 | const G4double ztmax = 0.99 ;
|
|---|
| 695 | G4double zt = zmean/tPathLength ;
|
|---|
| 696 |
|
|---|
| 697 | if (tPathLength > stepmin && zt < ztmax)
|
|---|
| 698 | {
|
|---|
| 699 | G4double u,cz1;
|
|---|
| 700 | if(zt >= third)
|
|---|
| 701 | {
|
|---|
| 702 | G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
|
|---|
| 703 | cz1 = 1.+cz ;
|
|---|
| 704 | G4double u0 = cz/cz1 ;
|
|---|
| 705 | G4double grej ;
|
|---|
| 706 | do {
|
|---|
| 707 | u = exp(log(G4UniformRand())/cz1) ;
|
|---|
| 708 | grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
|
|---|
| 709 | } while (grej < G4UniformRand()) ;
|
|---|
| 710 | }
|
|---|
| 711 | else
|
|---|
| 712 | {
|
|---|
| 713 | cz1 = 1./zt-1.;
|
|---|
| 714 | u = 1.-exp(log(G4UniformRand())/cz1) ;
|
|---|
| 715 | }
|
|---|
| 716 | zPathLength = tPathLength*u ;
|
|---|
| 717 | }
|
|---|
| 718 | }
|
|---|
| 719 |
|
|---|
| 720 | if(zPathLength > lambda0) zPathLength = lambda0;
|
|---|
| 721 | //G4cout << "zPathLength= " << zPathLength << " lambda1= " << lambda0 << G4endl;
|
|---|
| 722 | return zPathLength;
|
|---|
| 723 | }
|
|---|
| 724 |
|
|---|
| 725 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 726 |
|
|---|
| 727 | G4double G4UrbanMscModel93::ComputeTrueStepLength(G4double geomStepLength)
|
|---|
| 728 | {
|
|---|
| 729 | // step defined other than transportation
|
|---|
| 730 | if(geomStepLength == zPathLength && tPathLength <= currentRange)
|
|---|
| 731 | return tPathLength;
|
|---|
| 732 |
|
|---|
| 733 | // t = z for very small step
|
|---|
| 734 | zPathLength = geomStepLength;
|
|---|
| 735 | tPathLength = geomStepLength;
|
|---|
| 736 | if(geomStepLength < tlimitminfix) return tPathLength;
|
|---|
| 737 |
|
|---|
| 738 | // recalculation
|
|---|
| 739 | if((geomStepLength > lambda0*tausmall) && !insideskin)
|
|---|
| 740 | {
|
|---|
| 741 | if(par1 < 0.)
|
|---|
| 742 | tPathLength = -lambda0*log(1.-geomStepLength/lambda0) ;
|
|---|
| 743 | else
|
|---|
| 744 | {
|
|---|
| 745 | if(par1*par3*geomStepLength < 1.)
|
|---|
| 746 | tPathLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
|
|---|
| 747 | else
|
|---|
| 748 | tPathLength = currentRange;
|
|---|
| 749 | }
|
|---|
| 750 | }
|
|---|
| 751 | if(tPathLength < geomStepLength) tPathLength = geomStepLength;
|
|---|
| 752 | //G4cout << "tPathLength= " << tPathLength << " step= " << geomStepLength << G4endl;
|
|---|
| 753 |
|
|---|
| 754 | return tPathLength;
|
|---|
| 755 | }
|
|---|
| 756 |
|
|---|
| 757 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 758 |
|
|---|
| 759 | G4double G4UrbanMscModel93::ComputeTheta0(G4double trueStepLength,
|
|---|
| 760 | G4double KineticEnergy)
|
|---|
| 761 | {
|
|---|
| 762 | // for all particles take the width of the central part
|
|---|
| 763 | // from a parametrization similar to the Highland formula
|
|---|
| 764 | // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
|
|---|
| 765 | const G4double c_highland = 13.6*MeV ;
|
|---|
| 766 | G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
|
|---|
| 767 | KineticEnergy*(KineticEnergy+2.*mass)/
|
|---|
| 768 | ((currentKinEnergy+mass)*(KineticEnergy+mass)));
|
|---|
| 769 | y = trueStepLength/currentRadLength;
|
|---|
| 770 | G4double theta0 = c_highland*std::abs(charge)*sqrt(y)/betacp;
|
|---|
| 771 | y = log(y);
|
|---|
| 772 | // correction factor from e- scattering data
|
|---|
| 773 | G4double corr = coeffth1+coeffth2*y;
|
|---|
| 774 |
|
|---|
| 775 | theta0 *= corr ;
|
|---|
| 776 |
|
|---|
| 777 | return theta0;
|
|---|
| 778 | }
|
|---|
| 779 |
|
|---|
| 780 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 781 |
|
|---|
| 782 | void G4UrbanMscModel93::SampleScattering(const G4DynamicParticle* dynParticle,
|
|---|
| 783 | G4double safety)
|
|---|
| 784 | {
|
|---|
| 785 | G4double kineticEnergy = dynParticle->GetKineticEnergy();
|
|---|
| 786 |
|
|---|
| 787 | if((kineticEnergy <= 0.0) || (tPathLength <= tlimitminfix) ||
|
|---|
| 788 | (tPathLength/tausmall < lambda0)) return;
|
|---|
| 789 |
|
|---|
| 790 | G4double cth = SampleCosineTheta(tPathLength,kineticEnergy);
|
|---|
| 791 |
|
|---|
| 792 | // protection against 'bad' cth values
|
|---|
| 793 | if(std::abs(cth) > 1.) return;
|
|---|
| 794 |
|
|---|
| 795 | G4double sth = sqrt((1.0 - cth)*(1.0 + cth));
|
|---|
| 796 | G4double phi = twopi*G4UniformRand();
|
|---|
| 797 | G4double dirx = sth*cos(phi);
|
|---|
| 798 | G4double diry = sth*sin(phi);
|
|---|
| 799 |
|
|---|
| 800 | G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
|
|---|
| 801 | G4ThreeVector newDirection(dirx,diry,cth);
|
|---|
| 802 | newDirection.rotateUz(oldDirection);
|
|---|
| 803 | fParticleChange->ProposeMomentumDirection(newDirection);
|
|---|
| 804 |
|
|---|
| 805 | if (latDisplasment && safety > tlimitminfix) {
|
|---|
| 806 |
|
|---|
| 807 | G4double r = SampleDisplacement();
|
|---|
| 808 | /*
|
|---|
| 809 | G4cout << "G4UrbanMscModel93::SampleSecondaries: e(MeV)= " << kineticEnergy
|
|---|
| 810 | << " sinTheta= " << sth << " r(mm)= " << r
|
|---|
| 811 | << " trueStep(mm)= " << tPathLength
|
|---|
| 812 | << " geomStep(mm)= " << zPathLength
|
|---|
| 813 | << G4endl;
|
|---|
| 814 | */
|
|---|
| 815 | if(r > 0.)
|
|---|
| 816 | {
|
|---|
| 817 | G4double latcorr = LatCorrelation();
|
|---|
| 818 | if(latcorr > r) latcorr = r;
|
|---|
| 819 |
|
|---|
| 820 | // sample direction of lateral displacement
|
|---|
| 821 | // compute it from the lateral correlation
|
|---|
| 822 | G4double Phi = 0.;
|
|---|
| 823 | if(std::abs(r*sth) < latcorr)
|
|---|
| 824 | Phi = twopi*G4UniformRand();
|
|---|
| 825 | else
|
|---|
| 826 | {
|
|---|
| 827 | G4double psi = std::acos(latcorr/(r*sth));
|
|---|
| 828 | if(G4UniformRand() < 0.5)
|
|---|
| 829 | Phi = phi+psi;
|
|---|
| 830 | else
|
|---|
| 831 | Phi = phi-psi;
|
|---|
| 832 | }
|
|---|
| 833 |
|
|---|
| 834 | dirx = std::cos(Phi);
|
|---|
| 835 | diry = std::sin(Phi);
|
|---|
| 836 |
|
|---|
| 837 | G4ThreeVector latDirection(dirx,diry,0.0);
|
|---|
| 838 | latDirection.rotateUz(oldDirection);
|
|---|
| 839 |
|
|---|
| 840 | ComputeDisplacement(fParticleChange, latDirection, r, safety);
|
|---|
| 841 | }
|
|---|
| 842 | }
|
|---|
| 843 | }
|
|---|
| 844 |
|
|---|
| 845 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 846 |
|
|---|
| 847 | G4double G4UrbanMscModel93::SampleCosineTheta(G4double trueStepLength,
|
|---|
| 848 | G4double KineticEnergy)
|
|---|
| 849 | {
|
|---|
| 850 | G4double cth = 1. ;
|
|---|
| 851 | G4double tau = trueStepLength/lambda0 ;
|
|---|
| 852 |
|
|---|
| 853 | Zeff = couple->GetMaterial()->GetTotNbOfElectPerVolume()/
|
|---|
| 854 | couple->GetMaterial()->GetTotNbOfAtomsPerVolume() ;
|
|---|
| 855 |
|
|---|
| 856 | if(Zold != Zeff)
|
|---|
| 857 | UpdateCache();
|
|---|
| 858 |
|
|---|
| 859 | if(insideskin)
|
|---|
| 860 | {
|
|---|
| 861 | //no scattering, single or plural scattering
|
|---|
| 862 | G4double mean = trueStepLength/stepmin ;
|
|---|
| 863 |
|
|---|
| 864 | G4int n = G4Poisson(mean);
|
|---|
| 865 | if(n > 0)
|
|---|
| 866 | {
|
|---|
| 867 | //screening (Moliere-Bethe)
|
|---|
| 868 | G4double mom2 = KineticEnergy*(2.*mass+KineticEnergy);
|
|---|
| 869 | G4double beta2 = mom2/((KineticEnergy+mass)*(KineticEnergy+mass));
|
|---|
| 870 | G4double ascr = scr1/mom2;
|
|---|
| 871 | ascr *= 1.13+scr2/beta2;
|
|---|
| 872 | G4double ascr1 = 1.+2.*ascr;
|
|---|
| 873 | G4double bp1=ascr1+1.;
|
|---|
| 874 | G4double bm1=ascr1-1.;
|
|---|
| 875 |
|
|---|
| 876 | // single scattering from screened Rutherford x-section
|
|---|
| 877 | G4double ct,st,phi;
|
|---|
| 878 | G4double sx=0.,sy=0.,sz=0.;
|
|---|
| 879 | for(G4int i=1; i<=n; i++)
|
|---|
| 880 | {
|
|---|
| 881 | ct = ascr1-bp1*bm1/(2.*G4UniformRand()+bm1);
|
|---|
| 882 | if(ct < -1.) ct = -1.;
|
|---|
| 883 | if(ct > 1.) ct = 1.;
|
|---|
| 884 | st = sqrt(1.-ct*ct);
|
|---|
| 885 | phi = twopi*G4UniformRand();
|
|---|
| 886 | sx += st*cos(phi);
|
|---|
| 887 | sy += st*sin(phi);
|
|---|
| 888 | sz += ct;
|
|---|
| 889 | }
|
|---|
| 890 | cth = sz/sqrt(sx*sx+sy*sy+sz*sz);
|
|---|
| 891 | }
|
|---|
| 892 | }
|
|---|
| 893 | else
|
|---|
| 894 | {
|
|---|
| 895 | if(trueStepLength >= currentRange*dtrl)
|
|---|
| 896 | {
|
|---|
| 897 | if(par1*trueStepLength < 1.)
|
|---|
| 898 | tau = -par2*log(1.-par1*trueStepLength) ;
|
|---|
| 899 | // for the case if ioni/brems are inactivated
|
|---|
| 900 | // see the corresponding condition in ComputeGeomPathLength
|
|---|
| 901 | else if(1.-KineticEnergy/currentKinEnergy > taulim)
|
|---|
| 902 | tau = taubig ;
|
|---|
| 903 | }
|
|---|
| 904 | currentTau = tau ;
|
|---|
| 905 | lambdaeff = trueStepLength/currentTau;
|
|---|
| 906 | currentRadLength = couple->GetMaterial()->GetRadlen();
|
|---|
| 907 |
|
|---|
| 908 | if (tau >= taubig) cth = -1.+2.*G4UniformRand();
|
|---|
| 909 | else if (tau >= tausmall)
|
|---|
| 910 | {
|
|---|
| 911 | G4double xsi = 3.;
|
|---|
| 912 | G4double x0 = 1.;
|
|---|
| 913 | G4double a = 1., ea = 0., eaa = 1.;
|
|---|
| 914 | G4double b=2.,b1=3.,bx=1.,eb1=3.,ebx=1.;
|
|---|
| 915 | G4double prob = 1. , qprob = 1. ;
|
|---|
| 916 | G4double xmean1 = 1., xmean2 = 0.;
|
|---|
| 917 | G4double xmeanth = exp(-tau);
|
|---|
| 918 | G4double x2meanth = (1.+2.*exp(-2.5*tau))/3.;
|
|---|
| 919 |
|
|---|
| 920 | G4double relloss = 1.-KineticEnergy/currentKinEnergy;
|
|---|
| 921 | if(relloss > rellossmax)
|
|---|
| 922 | return SimpleScattering(xmeanth,x2meanth);
|
|---|
| 923 |
|
|---|
| 924 | G4double theta0 = ComputeTheta0(trueStepLength,KineticEnergy);
|
|---|
| 925 |
|
|---|
| 926 | // protection for very small angles
|
|---|
| 927 | if(theta0 < tausmall) return cth;
|
|---|
| 928 |
|
|---|
| 929 | if(theta0 > theta0max)
|
|---|
| 930 | return SimpleScattering(xmeanth,x2meanth);
|
|---|
| 931 | G4double sth = sin(0.5*theta0);
|
|---|
| 932 | a = 0.25/(sth*sth);
|
|---|
| 933 |
|
|---|
| 934 | ea = exp(-xsi);
|
|---|
| 935 | eaa = 1.-ea ;
|
|---|
| 936 | xmean1 = 1.-(1.-(1.+xsi)*ea)/(a*eaa);
|
|---|
| 937 | x0 = 1.-xsi/a;
|
|---|
| 938 |
|
|---|
| 939 | if(xmean1 <= 0.999*xmeanth)
|
|---|
| 940 | return SimpleScattering(xmeanth,x2meanth);
|
|---|
| 941 |
|
|---|
| 942 | // from e- and muon scattering data
|
|---|
| 943 | G4double c = coeffc1+coeffc2*y; ;
|
|---|
| 944 |
|
|---|
| 945 | if(abs(c-3.) < 0.001) c = 3.001;
|
|---|
| 946 | if(abs(c-2.) < 0.001) c = 2.001;
|
|---|
| 947 | if(abs(c-1.) < 0.001) c = 1.001;
|
|---|
| 948 |
|
|---|
| 949 | G4double c1 = c-1.;
|
|---|
| 950 |
|
|---|
| 951 | //from continuity of derivatives
|
|---|
| 952 | b = 1.+(c-xsi)/a;
|
|---|
| 953 |
|
|---|
| 954 | b1 = b+1.;
|
|---|
| 955 | bx = c/a;
|
|---|
| 956 | eb1 = exp(c1*log(b1));
|
|---|
| 957 | ebx = exp(c1*log(bx));
|
|---|
| 958 |
|
|---|
| 959 | xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/(c-2.))/(eb1-ebx);
|
|---|
| 960 |
|
|---|
| 961 | G4double f1x0 = a*ea/eaa;
|
|---|
| 962 | G4double f2x0 = c1*eb1/(bx*(eb1-ebx));
|
|---|
| 963 | prob = f2x0/(f1x0+f2x0);
|
|---|
| 964 |
|
|---|
| 965 | qprob = xmeanth/(prob*xmean1+(1.-prob)*xmean2);
|
|---|
| 966 |
|
|---|
| 967 | // sampling of costheta
|
|---|
| 968 | if(G4UniformRand() < qprob)
|
|---|
| 969 | {
|
|---|
| 970 | if(G4UniformRand() < prob)
|
|---|
| 971 | cth = 1.+log(ea+G4UniformRand()*eaa)/a ;
|
|---|
| 972 | else
|
|---|
| 973 | cth = b-b1*bx/exp(log(ebx+(eb1-ebx)*G4UniformRand())/c1) ;
|
|---|
| 974 | }
|
|---|
| 975 | else
|
|---|
| 976 | cth = -1.+2.*G4UniformRand();
|
|---|
| 977 | }
|
|---|
| 978 | }
|
|---|
| 979 | return cth ;
|
|---|
| 980 | }
|
|---|
| 981 |
|
|---|
| 982 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 983 |
|
|---|
| 984 | G4double G4UrbanMscModel93::SimpleScattering(G4double xmeanth,G4double x2meanth)
|
|---|
| 985 | {
|
|---|
| 986 | // 'large angle scattering'
|
|---|
| 987 | // 2 model functions with correct xmean and x2mean
|
|---|
| 988 | G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
|
|---|
| 989 | G4double prob = (a+2.)*xmeanth/a;
|
|---|
| 990 |
|
|---|
| 991 | // sampling
|
|---|
| 992 | G4double cth = 1.;
|
|---|
| 993 | if(G4UniformRand() < prob)
|
|---|
| 994 | cth = -1.+2.*exp(log(G4UniformRand())/(a+1.));
|
|---|
| 995 | else
|
|---|
| 996 | cth = -1.+2.*G4UniformRand();
|
|---|
| 997 | return cth;
|
|---|
| 998 | }
|
|---|
| 999 |
|
|---|
| 1000 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 1001 |
|
|---|
| 1002 | G4double G4UrbanMscModel93::SampleDisplacement()
|
|---|
| 1003 | {
|
|---|
| 1004 | const G4double kappa = 2.5;
|
|---|
| 1005 | const G4double kappapl1 = kappa+1.;
|
|---|
| 1006 | const G4double kappami1 = kappa-1.;
|
|---|
| 1007 | G4double rmean = 0.0;
|
|---|
| 1008 | if ((currentTau >= tausmall) && !insideskin) {
|
|---|
| 1009 | if (currentTau < taulim) {
|
|---|
| 1010 | rmean = kappa*currentTau*currentTau*currentTau*
|
|---|
| 1011 | (1.-kappapl1*currentTau*0.25)/6. ;
|
|---|
| 1012 |
|
|---|
| 1013 | } else {
|
|---|
| 1014 | G4double etau = 0.0;
|
|---|
| 1015 | if (currentTau<taubig) etau = exp(-currentTau);
|
|---|
| 1016 | rmean = -kappa*currentTau;
|
|---|
| 1017 | rmean = -exp(rmean)/(kappa*kappami1);
|
|---|
| 1018 | rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
|
|---|
| 1019 | }
|
|---|
| 1020 | if (rmean>0.) rmean = 2.*lambdaeff*sqrt(rmean/3.0);
|
|---|
| 1021 | else rmean = 0.;
|
|---|
| 1022 | }
|
|---|
| 1023 |
|
|---|
| 1024 | // protection against z > t ...........................
|
|---|
| 1025 | if(rmean > 0.) {
|
|---|
| 1026 | G4double zt = (tPathLength-zPathLength)*(tPathLength+zPathLength);
|
|---|
| 1027 | if(zt <= 0.)
|
|---|
| 1028 | rmean = 0.;
|
|---|
| 1029 | else if(rmean*rmean > zt)
|
|---|
| 1030 | rmean = sqrt(zt);
|
|---|
| 1031 | }
|
|---|
| 1032 | return rmean;
|
|---|
| 1033 | }
|
|---|
| 1034 |
|
|---|
| 1035 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 1036 |
|
|---|
| 1037 | G4double G4UrbanMscModel93::LatCorrelation()
|
|---|
| 1038 | {
|
|---|
| 1039 | const G4double kappa = 2.5;
|
|---|
| 1040 | const G4double kappami1 = kappa-1.;
|
|---|
| 1041 |
|
|---|
| 1042 | G4double latcorr = 0.;
|
|---|
| 1043 | if((currentTau >= tausmall) && !insideskin)
|
|---|
| 1044 | {
|
|---|
| 1045 | if(currentTau < taulim)
|
|---|
| 1046 | latcorr = lambdaeff*kappa*currentTau*currentTau*
|
|---|
| 1047 | (1.-(kappa+1.)*currentTau/3.)/3.;
|
|---|
| 1048 | else
|
|---|
| 1049 | {
|
|---|
| 1050 | G4double etau = 0.;
|
|---|
| 1051 | if(currentTau < taubig) etau = exp(-currentTau);
|
|---|
| 1052 | latcorr = -kappa*currentTau;
|
|---|
| 1053 | latcorr = exp(latcorr)/kappami1;
|
|---|
| 1054 | latcorr += 1.-kappa*etau/kappami1 ;
|
|---|
| 1055 | latcorr *= 2.*lambdaeff/3. ;
|
|---|
| 1056 | }
|
|---|
| 1057 | }
|
|---|
| 1058 |
|
|---|
| 1059 | return latcorr;
|
|---|
| 1060 | }
|
|---|
| 1061 |
|
|---|
| 1062 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|