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