[819] | 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 | // |
---|
[1196] | 26 | // $Id: G4MscModel71.cc,v 1.8 2009/11/01 13:05:01 vnivanch Exp $ |
---|
[1228] | 27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4MscModel71 |
---|
| 35 | // |
---|
| 36 | // Author: Laszlo Urban |
---|
| 37 | // |
---|
| 38 | // Creation date: 03.03.2001 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 27-03-03 Move model part from G4MultipleScattering (V.Ivanchenko) |
---|
| 43 | // 23-05-03 important change in angle distribution for muons/hadrons |
---|
| 44 | // the central part now is similar to the Highland parametrization + |
---|
| 45 | // minor correction in angle sampling algorithm (for all particles) |
---|
| 46 | // (L.Urban) |
---|
| 47 | // 30-05-03 misprint in SampleCosineTheta corrected(L.Urban) |
---|
| 48 | // 27-03-03 Rename (V.Ivanchenko) |
---|
| 49 | // 05-08-03 angle distribution has been modified (L.Urban) |
---|
| 50 | // 06-11-03 precision problems solved for high energy (PeV) particles |
---|
| 51 | // change in the tail of the angular distribution |
---|
| 52 | // highKinEnergy is set to 100 PeV (L.Urban) |
---|
| 53 | // |
---|
| 54 | // 10-11-03 highKinEnergy is set back to 100 TeV, some tail tuning + |
---|
| 55 | // cleaning (L.Urban) |
---|
| 56 | // 26-11-03 correction in TrueStepLength : |
---|
| 57 | // trueLength <= currentRange (L.Urban) |
---|
| 58 | // 01-03-04 signature changed in SampleCosineTheta, |
---|
| 59 | // energy dependence calculations has been simplified, |
---|
| 60 | // 11-03-04 corrections in GeomPathLength,TrueStepLength, |
---|
| 61 | // SampleCosineTheta |
---|
| 62 | // 23-04-04 true -> geom and geom -> true transformation has been |
---|
| 63 | // rewritten, changes in the angular distribution (L.Urban) |
---|
| 64 | // 19-07-04 correction in SampleCosineTheta in order to avoid |
---|
| 65 | // num. precision problems at high energy/small step(L.Urban) |
---|
| 66 | // 17-08-04 changes in the angle distribution (slightly modified |
---|
| 67 | // Highland formula for the width of the central part, |
---|
| 68 | // changes in the numerical values of some other parameters) |
---|
| 69 | // ---> approximately step independent distribution (L.Urban) |
---|
| 70 | // 21-09-04 change in the tail of the angular distribution (L.Urban) |
---|
| 71 | // |
---|
| 72 | // 03-11-04 precision problem for very high energy ions and small stepsize |
---|
| 73 | // solved in SampleCosineTheta (L.Urban). |
---|
| 74 | // 15-04-05 optimize internal interface - add SampleSecondaries method (V.Ivanchenko) |
---|
| 75 | // 03-10-05 Model is freezed with the name McsModel71 (V.Ivanchenko) |
---|
| 76 | // 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko) |
---|
| 77 | // |
---|
| 78 | |
---|
| 79 | // Class Description: |
---|
| 80 | // |
---|
| 81 | // Implementation of the model of multiple scattering based on |
---|
| 82 | // H.W.Lewis Phys Rev 78 (1950) 526 and others |
---|
| 83 | |
---|
| 84 | // ------------------------------------------------------------------- |
---|
| 85 | // |
---|
| 86 | |
---|
| 87 | |
---|
| 88 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 89 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 90 | |
---|
| 91 | #include "G4MscModel71.hh" |
---|
| 92 | #include "Randomize.hh" |
---|
| 93 | #include "G4Electron.hh" |
---|
| 94 | #include "G4LossTableManager.hh" |
---|
| 95 | #include "G4PhysicsTable.hh" |
---|
| 96 | #include "G4ParticleChangeForMSC.hh" |
---|
| 97 | #include "G4TransportationManager.hh" |
---|
| 98 | #include "G4Navigator.hh" |
---|
| 99 | |
---|
| 100 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 101 | |
---|
| 102 | using namespace std; |
---|
| 103 | |
---|
| 104 | G4MscModel71::G4MscModel71(G4double& m_dtrl, G4double& m_NuclCorrPar, |
---|
| 105 | G4double& m_FactPar, G4double& m_factail, |
---|
[1196] | 106 | G4bool& m_samplez, const G4String& nam) |
---|
[819] | 107 | : G4VEmModel(nam), |
---|
| 108 | taubig(8.0), |
---|
| 109 | tausmall(1.e-20), |
---|
| 110 | taulim(1.e-6), |
---|
| 111 | dtrl(m_dtrl), |
---|
| 112 | NuclCorrPar (m_NuclCorrPar), |
---|
| 113 | FactPar(m_FactPar), |
---|
| 114 | factail(m_factail), |
---|
| 115 | samplez(m_samplez), |
---|
| 116 | isInitialized(false) |
---|
| 117 | { |
---|
| 118 | stepmin = 1.e-6*mm; |
---|
[1196] | 119 | currentRange = 0.; |
---|
| 120 | G4cout << G4endl; |
---|
| 121 | G4cout << "!!! G4MscModel71 class is obsolete and will be removed for the next major Geant4 release !!!" << G4endl; |
---|
| 122 | G4cout << "!!! Please use other models (G4UrbanMscModel90, 92, 93) !!!" << G4endl; |
---|
| 123 | G4cout << G4endl; |
---|
[819] | 124 | } |
---|
| 125 | |
---|
| 126 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 127 | |
---|
| 128 | G4MscModel71::~G4MscModel71() |
---|
| 129 | {} |
---|
| 130 | |
---|
| 131 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 132 | |
---|
| 133 | void G4MscModel71::Initialise(const G4ParticleDefinition* p, |
---|
| 134 | const G4DataVector&) |
---|
| 135 | { |
---|
| 136 | if(isInitialized) return; |
---|
| 137 | // set values of some data members |
---|
| 138 | sigmafactor = twopi*classic_electr_radius*classic_electr_radius; |
---|
| 139 | particle = p; |
---|
| 140 | mass = particle->GetPDGMass(); |
---|
| 141 | charge = particle->GetPDGCharge()/eplus; |
---|
| 142 | b = 1. ; |
---|
| 143 | xsi = 3.00 ; |
---|
| 144 | |
---|
| 145 | if(pParticleChange) |
---|
[1055] | 146 | fParticleChange = static_cast<G4ParticleChangeForMSC*>(pParticleChange); |
---|
[819] | 147 | else |
---|
| 148 | fParticleChange = new G4ParticleChangeForMSC(); |
---|
| 149 | |
---|
| 150 | navigator = G4TransportationManager::GetTransportationManager() |
---|
| 151 | ->GetNavigatorForTracking(); |
---|
| 152 | } |
---|
| 153 | |
---|
| 154 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 155 | |
---|
| 156 | G4double G4MscModel71::ComputeCrossSectionPerAtom( |
---|
| 157 | const G4ParticleDefinition* part, |
---|
| 158 | G4double KineticEnergy, |
---|
| 159 | G4double AtomicNumber, |
---|
| 160 | G4double AtomicWeight, |
---|
| 161 | G4double, |
---|
| 162 | G4double) |
---|
| 163 | { |
---|
| 164 | const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2* |
---|
| 165 | Bohr_radius*Bohr_radius/(hbarc*hbarc); |
---|
| 166 | const G4double epsmin = 1.e-4 , epsmax = 1.e10; |
---|
| 167 | |
---|
| 168 | const G4double Zdat[15] = { 4., 6.,13.,20.,26.,29.,32.,38.,47., |
---|
| 169 | 50.,56.,64.,74.,79.,82. }; |
---|
| 170 | |
---|
| 171 | const G4double Tdat[23] = {0.0001*MeV,0.0002*MeV,0.0004*MeV,0.0007*MeV, |
---|
| 172 | 0.001*MeV,0.002*MeV,0.004*MeV,0.007*MeV, |
---|
| 173 | 0.01*MeV,0.02*MeV,0.04*MeV,0.07*MeV, |
---|
| 174 | 0.1*MeV,0.2*MeV,0.4*MeV,0.7*MeV, |
---|
| 175 | 1.*MeV,2.*MeV,4.*MeV,7.*MeV,10.*MeV,20.*MeV, |
---|
| 176 | 10000.0*MeV}; |
---|
| 177 | |
---|
| 178 | // corr. factors for e-/e+ lambda |
---|
| 179 | G4double celectron[15][23] = |
---|
| 180 | {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054, |
---|
| 181 | 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111, |
---|
| 182 | 1.112,1.108,1.100,1.093,1.089,1.087,0.7235 }, |
---|
| 183 | {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051, |
---|
| 184 | 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108, |
---|
| 185 | 1.109,1.105,1.097,1.090,1.086,1.082,0.7925 }, |
---|
| 186 | {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156, |
---|
| 187 | 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132, |
---|
| 188 | 1.131,1.124,1.113,1.104,1.099,1.098,0.9147 }, |
---|
| 189 | {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236, |
---|
| 190 | 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113, |
---|
| 191 | 1.112,1.105,1.096,1.089,1.085,1.098,0.9700 }, |
---|
| 192 | {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265, |
---|
| 193 | 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073, |
---|
| 194 | 1.073,1.070,1.064,1.059,1.056,1.056,1.0022 }, |
---|
| 195 | {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330, |
---|
| 196 | 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074, |
---|
| 197 | 1.074,1.070,1.063,1.059,1.056,1.052,1.0158 }, |
---|
| 198 | {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386, |
---|
| 199 | 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069, |
---|
| 200 | 1.068,1.064,1.059,1.054,1.051,1.050,1.0284 }, |
---|
| 201 | {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439, |
---|
| 202 | 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039, |
---|
| 203 | 1.039,1.037,1.034,1.031,1.030,1.036,1.0515 }, |
---|
| 204 | {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631, |
---|
| 205 | 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033, |
---|
| 206 | 1.031,1.028,1.024,1.022,1.021,1.024,1.0834 }, |
---|
| 207 | {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669, |
---|
| 208 | 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022, |
---|
| 209 | 1.020,1.017,1.015,1.013,1.013,1.020,1.0937 }, |
---|
| 210 | {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720, |
---|
| 211 | 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997, |
---|
| 212 | 0.995,0.993,0.993,0.993,0.993,1.011,1.1140 }, |
---|
| 213 | {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855, |
---|
| 214 | 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976, |
---|
| 215 | 0.974,0.972,0.973,0.974,0.975,0.987,1.1410 }, |
---|
| 216 | {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059, |
---|
| 217 | 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954, |
---|
| 218 | 0.950,0.947,0.949,0.952,0.954,0.963,1.1750 }, |
---|
| 219 | {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182, |
---|
| 220 | 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947, |
---|
| 221 | 0.941,0.938,0.940,0.944,0.946,0.954,1.1922 }, |
---|
| 222 | {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239, |
---|
| 223 | 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939, |
---|
| 224 | 0.933,0.930,0.933,0.936,0.939,0.949,1.2026 }}; |
---|
| 225 | G4double cpositron[15][23] = { |
---|
| 226 | {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110, |
---|
| 227 | 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131, |
---|
| 228 | 1.131,1.126,1.117,1.108,1.103,1.100,0.7235 }, |
---|
| 229 | {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145, |
---|
| 230 | 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137, |
---|
| 231 | 1.138,1.132,1.122,1.113,1.108,1.102,0.7925 }, |
---|
| 232 | {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451, |
---|
| 233 | 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205, |
---|
| 234 | 1.203,1.190,1.173,1.159,1.151,1.145,0.9147 }, |
---|
| 235 | {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715, |
---|
| 236 | 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228, |
---|
| 237 | 1.225,1.210,1.191,1.175,1.166,1.174,0.9700 }, |
---|
| 238 | {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820, |
---|
| 239 | 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219, |
---|
| 240 | 1.217,1.203,1.184,1.169,1.160,1.151,1.0022 }, |
---|
| 241 | {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996, |
---|
| 242 | 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241, |
---|
| 243 | 1.237,1.222,1.201,1.184,1.174,1.159,1.0158 }, |
---|
| 244 | {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155, |
---|
| 245 | 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256, |
---|
| 246 | 1.252,1.234,1.212,1.194,1.183,1.170,1.0284 }, |
---|
| 247 | {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348, |
---|
| 248 | 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258, |
---|
| 249 | 1.254,1.237,1.214,1.195,1.185,1.179,1.0515 }, |
---|
| 250 | {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808, |
---|
| 251 | 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320, |
---|
| 252 | 1.312,1.288,1.258,1.235,1.221,1.205,1.0834 }, |
---|
| 253 | {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917, |
---|
| 254 | 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327, |
---|
| 255 | 1.320,1.294,1.264,1.240,1.226,1.214,1.0937 }, |
---|
| 256 | {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066, |
---|
| 257 | 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336, |
---|
| 258 | 1.328,1.302,1.270,1.245,1.231,1.233,1.1140 }, |
---|
| 259 | {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498, |
---|
| 260 | 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371, |
---|
| 261 | 1.361,1.330,1.294,1.267,1.251,1.239,1.1410 }, |
---|
| 262 | {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155, |
---|
| 263 | 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423, |
---|
| 264 | 1.409,1.372,1.330,1.298,1.280,1.258,1.1750 }, |
---|
| 265 | {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407, |
---|
| 266 | 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459, |
---|
| 267 | 1.442,1.400,1.354,1.319,1.299,1.272,1.1922 }, |
---|
| 268 | {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542, |
---|
| 269 | 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474, |
---|
| 270 | 1.456,1.412,1.364,1.328,1.307,1.282,1.2026 }}; |
---|
| 271 | |
---|
| 272 | G4double sigma; |
---|
| 273 | if (part != particle ) { |
---|
| 274 | particle = part; |
---|
| 275 | mass = particle->GetPDGMass(); |
---|
| 276 | charge = particle->GetPDGCharge()/eplus; |
---|
| 277 | } |
---|
| 278 | |
---|
| 279 | G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23); |
---|
| 280 | |
---|
| 281 | // correction if particle .ne. e-/e+ |
---|
| 282 | // compute equivalent kinetic energy |
---|
| 283 | // lambda depends on p*beta .... |
---|
| 284 | |
---|
| 285 | G4double eKineticEnergy = KineticEnergy; |
---|
| 286 | |
---|
| 287 | if((particle->GetParticleName() != "e-") && |
---|
| 288 | (particle->GetParticleName() != "e+") ) |
---|
| 289 | { |
---|
| 290 | G4double TAU = KineticEnergy/mass ; |
---|
| 291 | G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ; |
---|
| 292 | G4double w = c-2. ; |
---|
| 293 | G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ; |
---|
| 294 | eKineticEnergy = electron_mass_c2*tau ; |
---|
| 295 | } |
---|
| 296 | |
---|
| 297 | G4double ChargeSquare = charge*charge; |
---|
| 298 | |
---|
| 299 | G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ; |
---|
| 300 | G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2) |
---|
| 301 | /(eTotalEnergy*eTotalEnergy); |
---|
| 302 | G4double bg2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2) |
---|
| 303 | /(electron_mass_c2*electron_mass_c2); |
---|
| 304 | |
---|
| 305 | G4double eps = epsfactor*bg2/Z23; |
---|
| 306 | |
---|
| 307 | if (eps<epsmin) sigma = 2.*eps*eps; |
---|
| 308 | else if(eps<epsmax) sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps); |
---|
| 309 | else sigma = log(2.*eps)-1.+1./eps; |
---|
| 310 | |
---|
| 311 | sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2); |
---|
| 312 | |
---|
| 313 | // nuclear size effect correction for high energy |
---|
| 314 | // ( a simple approximation at present) |
---|
| 315 | G4double corrnuclsize,a,w1,w2,w; |
---|
| 316 | |
---|
| 317 | G4double x0 = 1. - NuclCorrPar*mass/(KineticEnergy* |
---|
| 318 | exp(log(AtomicWeight/(g/mole))/3.)); |
---|
| 319 | if ( x0 < -1. || eKineticEnergy <= 10.*MeV) |
---|
| 320 | { |
---|
| 321 | x0 = -1.; |
---|
| 322 | corrnuclsize = 1.; |
---|
| 323 | } |
---|
| 324 | else |
---|
| 325 | { |
---|
| 326 | a = 1.+1./eps; |
---|
| 327 | if (eps > epsmax) w1=log(2.*eps)+1./eps-3./(8.*eps*eps); |
---|
| 328 | else w1=log((a+1.)/(a-1.))-2./(a+1.); |
---|
| 329 | w = 1./((1.-x0)*eps); |
---|
| 330 | if (w < epsmin) w2=-log(w)-1.+2.*w-1.5*w*w; |
---|
| 331 | else w2 = log((a-x0)/(a-1.))-(1.-x0)/(a-x0); |
---|
| 332 | corrnuclsize = w1/w2; |
---|
| 333 | corrnuclsize = exp(-FactPar*mass/KineticEnergy)* |
---|
| 334 | (corrnuclsize-1.)+1.; |
---|
| 335 | } |
---|
| 336 | |
---|
| 337 | // interpolate in AtomicNumber and beta2 |
---|
| 338 | // get bin number in Z |
---|
| 339 | G4int iZ = 14; |
---|
| 340 | while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1; |
---|
| 341 | if (iZ==14) iZ = 13; |
---|
| 342 | if (iZ==-1) iZ = 0 ; |
---|
| 343 | |
---|
| 344 | G4double Z1 = Zdat[iZ]; |
---|
| 345 | G4double Z2 = Zdat[iZ+1]; |
---|
| 346 | G4double ratZ = (AtomicNumber-Z1)/(Z2-Z1); |
---|
| 347 | |
---|
| 348 | // get bin number in T (beta2) |
---|
| 349 | G4int iT = 22; |
---|
| 350 | while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1; |
---|
| 351 | if(iT==22) iT = 21; |
---|
| 352 | if(iT==-1) iT = 0 ; |
---|
| 353 | |
---|
| 354 | // calculate betasquare values |
---|
| 355 | G4double T = Tdat[iT], E = T + electron_mass_c2; |
---|
| 356 | G4double b2small = T*(E+electron_mass_c2)/(E*E); |
---|
| 357 | T = Tdat[iT+1]; E = T + electron_mass_c2; |
---|
| 358 | G4double b2big = T*(E+electron_mass_c2)/(E*E); |
---|
| 359 | G4double ratb2 = (beta2-b2small)/(b2big-b2small); |
---|
| 360 | G4double c1,c2,cc1,cc2,corr; |
---|
| 361 | if (charge < 0.) |
---|
| 362 | { |
---|
| 363 | c1 = celectron[iZ][iT]; |
---|
| 364 | c2 = celectron[iZ+1][iT]; |
---|
| 365 | cc1 = c1+ratZ*(c2-c1); |
---|
| 366 | |
---|
| 367 | c1 = celectron[iZ][iT+1]; |
---|
| 368 | c2 = celectron[iZ+1][iT+1]; |
---|
| 369 | cc2 = c1+ratZ*(c2-c1); |
---|
| 370 | |
---|
| 371 | corr = cc1+ratb2*(cc2-cc1); |
---|
| 372 | sigma /= corr; |
---|
| 373 | } |
---|
| 374 | |
---|
| 375 | if (charge > 0.) |
---|
| 376 | { |
---|
| 377 | c1 = cpositron[iZ][iT]; |
---|
| 378 | c2 = cpositron[iZ+1][iT]; |
---|
| 379 | cc1 = c1+ratZ*(c2-c1); |
---|
| 380 | |
---|
| 381 | c1 = cpositron[iZ][iT+1]; |
---|
| 382 | c2 = cpositron[iZ+1][iT+1]; |
---|
| 383 | cc2 = c1+ratZ*(c2-c1); |
---|
| 384 | |
---|
| 385 | corr = cc1+ratb2*(cc2-cc1); |
---|
| 386 | sigma /= corr; |
---|
| 387 | } |
---|
| 388 | |
---|
| 389 | sigma *= sigmafactor; |
---|
| 390 | |
---|
| 391 | // nucl. size correction for particles other than e+/e- only at present !!!! |
---|
| 392 | if((particle->GetParticleName() != "e-") && |
---|
| 393 | (particle->GetParticleName() != "e+") ) |
---|
| 394 | sigma /= corrnuclsize; |
---|
| 395 | // G4cout << "e= " << KineticEnergy << " sigma= " << sigma << G4endl; |
---|
| 396 | |
---|
| 397 | return sigma; |
---|
| 398 | } |
---|
| 399 | |
---|
| 400 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 401 | |
---|
| 402 | G4double G4MscModel71::GeomPathLength( |
---|
| 403 | G4PhysicsTable* theLambdaTable, |
---|
| 404 | const G4MaterialCutsCouple* couple, |
---|
| 405 | const G4ParticleDefinition* theParticle, |
---|
| 406 | G4double& T0, |
---|
| 407 | G4double lambda, |
---|
| 408 | G4double range, |
---|
| 409 | G4double truePathLength) |
---|
| 410 | { |
---|
| 411 | // do the true -> geom transformation |
---|
| 412 | const G4double ztmax = 101./103. ; |
---|
| 413 | if (theParticle != particle ) { |
---|
| 414 | particle = theParticle; |
---|
| 415 | mass = particle->GetPDGMass(); |
---|
| 416 | charge = particle->GetPDGCharge()/eplus; |
---|
| 417 | } |
---|
| 418 | currentKinEnergy = T0; |
---|
| 419 | currentRange = range ; |
---|
| 420 | currentRadLength = couple->GetMaterial()->GetRadlen(); |
---|
| 421 | |
---|
| 422 | lambda0 = lambda; |
---|
| 423 | par1 = -1. ; |
---|
| 424 | par2 = par3 = 0. ; |
---|
| 425 | tPathLength = truePathLength; |
---|
| 426 | |
---|
| 427 | // this correction needed to run MSC with eIoni and eBrem inactivated |
---|
| 428 | // and makes no harm for a normal run |
---|
| 429 | if(tPathLength > range) |
---|
| 430 | tPathLength = range ; |
---|
| 431 | |
---|
| 432 | G4double tau = tPathLength/lambda0 ; |
---|
| 433 | |
---|
| 434 | if (tau <= tausmall) return tPathLength; |
---|
| 435 | |
---|
| 436 | G4double zmean = tPathLength; |
---|
| 437 | if (tPathLength < range*dtrl) { |
---|
| 438 | zmean = lambda0*(1.-exp(-tau)); |
---|
| 439 | if(tau < taulim) zmean = tPathLength*(1.-0.5*tPathLength/lambda0) ; |
---|
| 440 | } else if(T0 < mass) { |
---|
| 441 | par1 = 1./range ; |
---|
| 442 | par2 = 1./(par1*lambda0) ; |
---|
| 443 | par3 = 1.+par2 ; |
---|
| 444 | zmean = (1.-exp(par3*log(1.-tPathLength/range)))/(par1*par3) ; |
---|
| 445 | } else { |
---|
| 446 | G4LossTableManager* theManager = G4LossTableManager::Instance(); |
---|
| 447 | G4double T1 = theManager->GetEnergy(particle,range-tPathLength,couple); |
---|
| 448 | G4double lambda1 ; |
---|
| 449 | if (theLambdaTable) { |
---|
| 450 | G4bool bb; |
---|
| 451 | lambda1 = ((*theLambdaTable)[couple->GetIndex()])->GetValue(T1,bb); |
---|
| 452 | } else { |
---|
| 453 | lambda1 = CrossSection(couple,particle,T1,0.0,1.0); |
---|
| 454 | } |
---|
| 455 | lambda1 = 1.0/lambda1; |
---|
| 456 | par1 = (lambda0-lambda1)/(lambda0*tPathLength) ; |
---|
| 457 | par2 = 1./(par1*lambda0) ; |
---|
| 458 | par3 = 1.+par2 ; |
---|
| 459 | zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ; |
---|
| 460 | } |
---|
| 461 | |
---|
| 462 | // sample z |
---|
| 463 | G4double zPathLength = zmean ; |
---|
| 464 | G4double zt = zmean/tPathLength ; |
---|
| 465 | if (tPathLength >= stepmin && samplez && zt > 0.5 && zt < ztmax) |
---|
| 466 | { |
---|
| 467 | G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ; |
---|
| 468 | G4double cz1 = 1.+cz ; |
---|
| 469 | G4double u0 = cz/cz1 ; |
---|
| 470 | G4double u,grej ; |
---|
| 471 | do { |
---|
| 472 | u = exp(log(G4UniformRand())/cz1) ; |
---|
| 473 | grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ; |
---|
| 474 | } while (grej < G4UniformRand()) ; |
---|
| 475 | zPathLength = tPathLength*u ; |
---|
| 476 | } |
---|
| 477 | |
---|
| 478 | return zPathLength ; |
---|
| 479 | } |
---|
| 480 | |
---|
| 481 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 482 | |
---|
| 483 | G4double G4MscModel71::TrueStepLength(G4double geomStepLength) |
---|
| 484 | { |
---|
| 485 | G4double trueLength = geomStepLength; |
---|
| 486 | trueLength = geomStepLength; |
---|
| 487 | if(geomStepLength > lambda0*tausmall) |
---|
| 488 | { |
---|
| 489 | if(par1 < 0.) |
---|
| 490 | trueLength = -lambda0*log(1.-geomStepLength/lambda0) ; |
---|
| 491 | else |
---|
| 492 | { |
---|
| 493 | if(par1*par3*geomStepLength < 1.) |
---|
| 494 | trueLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ; |
---|
| 495 | else |
---|
| 496 | trueLength = currentRange ; |
---|
| 497 | } |
---|
| 498 | } |
---|
| 499 | |
---|
| 500 | if(trueLength > tPathLength) trueLength = tPathLength; |
---|
| 501 | if(trueLength > currentRange) trueLength = currentRange ; |
---|
| 502 | if(trueLength < geomStepLength) trueLength = geomStepLength; |
---|
| 503 | |
---|
| 504 | return trueLength; |
---|
| 505 | } |
---|
| 506 | |
---|
| 507 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 508 | |
---|
| 509 | void G4MscModel71::SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 510 | const G4MaterialCutsCouple*, |
---|
| 511 | const G4DynamicParticle* dynParticle, |
---|
| 512 | G4double truestep, |
---|
| 513 | G4double safety) |
---|
| 514 | { |
---|
| 515 | G4double kineticEnergy = dynParticle->GetKineticEnergy(); |
---|
| 516 | if(kineticEnergy <= 0.0) return; |
---|
| 517 | |
---|
| 518 | G4double cth = SampleCosineTheta(truestep,kineticEnergy); |
---|
| 519 | G4double sth = sqrt((1.0 - cth)*(1.0 + cth)); |
---|
| 520 | G4double phi = twopi*G4UniformRand(); |
---|
| 521 | G4double dirx = sth*cos(phi); |
---|
| 522 | G4double diry = sth*sin(phi); |
---|
| 523 | |
---|
| 524 | G4ThreeVector oldDirection = dynParticle->GetMomentumDirection(); |
---|
| 525 | G4ThreeVector newDirection(dirx,diry,cth); |
---|
| 526 | newDirection.rotateUz(oldDirection); |
---|
| 527 | fParticleChange->ProposeMomentumDirection(newDirection); |
---|
| 528 | |
---|
| 529 | /* |
---|
| 530 | const G4ParticleDefinition* pd = dynParticle->GetDefinition(); |
---|
| 531 | G4cout << "G4MscModel71: Sample secondary; E(MeV)= " << kineticEnergy/MeV |
---|
| 532 | << " MeV; step(mm)= " << truestep/mm |
---|
| 533 | << ", safety(mm)= " << safety/mm << " " << pd->GetParticleName() |
---|
| 534 | << G4endl; |
---|
| 535 | */ |
---|
| 536 | |
---|
| 537 | if (latDisplasment && safety > 0.0) { |
---|
| 538 | |
---|
| 539 | G4double r = SampleDisplacement(); |
---|
| 540 | if (r > safety) r = safety; |
---|
| 541 | |
---|
| 542 | // sample direction of lateral displacement |
---|
| 543 | G4double phi = twopi*G4UniformRand(); |
---|
| 544 | G4double dirx = std::cos(phi); |
---|
| 545 | G4double diry = std::sin(phi); |
---|
| 546 | |
---|
| 547 | G4ThreeVector newPosition(dirx,diry,0.0); |
---|
| 548 | newPosition.rotateUz(oldDirection); |
---|
| 549 | |
---|
| 550 | // compute new endpoint of the Step |
---|
| 551 | newPosition *= r; |
---|
| 552 | newPosition += *(fParticleChange->GetProposedPosition()); |
---|
| 553 | |
---|
| 554 | navigator->LocateGlobalPointWithinVolume(newPosition); |
---|
| 555 | |
---|
| 556 | fParticleChange->ProposePosition(newPosition); |
---|
| 557 | } |
---|
| 558 | } |
---|
| 559 | |
---|
| 560 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 561 | |
---|
| 562 | G4double G4MscModel71::SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy) |
---|
| 563 | { |
---|
| 564 | G4double cth = 1. ; |
---|
| 565 | G4double tau = trueStepLength/lambda0 ; |
---|
| 566 | |
---|
[961] | 567 | if(trueStepLength >= currentRange*dtrl) { |
---|
[819] | 568 | if(par1*trueStepLength < 1.) |
---|
| 569 | tau = -par2*log(1.-par1*trueStepLength) ; |
---|
| 570 | else |
---|
| 571 | tau = taubig ; |
---|
[961] | 572 | } |
---|
[819] | 573 | currentTau = tau ; |
---|
| 574 | |
---|
| 575 | if(trueStepLength < stepmin) |
---|
| 576 | cth = exp(-tau) ; |
---|
| 577 | else |
---|
| 578 | { |
---|
| 579 | if (tau >= taubig) cth = -1.+2.*G4UniformRand(); |
---|
| 580 | else if (tau >= tausmall) |
---|
| 581 | { |
---|
| 582 | G4double a ; |
---|
| 583 | |
---|
| 584 | // for all particles take the width of the central part |
---|
| 585 | // from a parametrization similar to the Highland formula |
---|
| 586 | // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10) |
---|
| 587 | // here : theta0 = 13.6*MeV*Q*(t/X0)**0.555/(beta*cp) |
---|
| 588 | const G4double c_highland = 13.6*MeV, corr_highland=0.555 ; |
---|
| 589 | G4double Q = std::abs(charge) ; |
---|
| 590 | G4double xx0 = trueStepLength/currentRadLength; |
---|
| 591 | G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)* |
---|
| 592 | KineticEnergy*(KineticEnergy+2.*mass)/ |
---|
| 593 | ((currentKinEnergy+mass)*(KineticEnergy+mass))) ; |
---|
| 594 | G4double theta0 = c_highland*Q*exp(corr_highland*log(xx0))/betacp ; |
---|
| 595 | |
---|
| 596 | if(theta0 > taulim) a = 0.5/(1.-cos(theta0)) ; |
---|
| 597 | else a = 1.0/(theta0*theta0) ; |
---|
| 598 | |
---|
| 599 | G4double xmeanth = exp(-tau); |
---|
| 600 | G4double xmeanth1 = 1.-xmeanth ; |
---|
| 601 | if(currentTau < taulim) xmeanth1 = tau ; |
---|
| 602 | |
---|
| 603 | const G4double x1fac1 = exp(-xsi) ; |
---|
| 604 | const G4double x1fac2 = (1.-(1.+xsi)*x1fac1)/(1.-x1fac1) ; |
---|
| 605 | const G4double x1fac3 = 1.3 ; |
---|
| 606 | |
---|
| 607 | G4double ea,eaa,xmean1 ; |
---|
| 608 | G4double c = 2.,b1 = 2., bx = 2., |
---|
| 609 | eb1 = b1, ebx = b1, xmean2 = 0. ; |
---|
| 610 | G4double prob = 1., qprob ; |
---|
| 611 | G4double x0 = 1.-xsi/a; |
---|
| 612 | G4double oneminusx0=xsi/a ; |
---|
| 613 | G4double oneplusx0=2.+xsi/a ; |
---|
| 614 | |
---|
| 615 | G4double f1x0=1., f2x0=1. ; |
---|
| 616 | const G4double tau0 = 0.10 ; |
---|
| 617 | if(tau > tau0) |
---|
| 618 | { |
---|
| 619 | // 1 model function |
---|
| 620 | a = 1./xmeanth1 ; |
---|
| 621 | ea = exp(-2.*a) ; |
---|
| 622 | eaa= 1.-ea ; |
---|
| 623 | xmean1 = 1.-1./a+2.*ea/eaa ; |
---|
| 624 | prob = 1. ; |
---|
| 625 | qprob = 1. ; |
---|
| 626 | } |
---|
| 627 | else if (x0 <= -1.) |
---|
| 628 | { |
---|
| 629 | // 2 model fuctions only |
---|
| 630 | // in order to have xmean1 > xmeanth -> qprob < 1 |
---|
| 631 | x0 = -1.; |
---|
| 632 | |
---|
| 633 | if( a < 1./xmeanth1) |
---|
| 634 | a = 1./xmeanth1 ; |
---|
| 635 | |
---|
| 636 | oneminusx0 = 1.-x0 ; |
---|
| 637 | oneplusx0 = 1.+x0 ; |
---|
| 638 | ea = exp(-a*oneminusx0); |
---|
| 639 | eaa = 1.-ea ; |
---|
| 640 | xmean1 = 1.-1./a+oneminusx0*ea/eaa ; |
---|
| 641 | qprob = xmeanth/xmean1 ; |
---|
| 642 | |
---|
| 643 | } |
---|
| 644 | else |
---|
| 645 | { |
---|
| 646 | // 3 model fuctions |
---|
| 647 | // in order to have xmean1 > xmeanth |
---|
| 648 | if((1.-x1fac2/a) < xmeanth) |
---|
| 649 | { |
---|
| 650 | a = x1fac3*x1fac2/xmeanth1 ; |
---|
| 651 | x0 = 1.-xsi/a ; |
---|
| 652 | oneminusx0=xsi/a ; |
---|
| 653 | oneplusx0=2.-xsi/a ; |
---|
| 654 | } |
---|
| 655 | |
---|
| 656 | ea = x1fac1 ; |
---|
| 657 | eaa = 1.-ea ; |
---|
| 658 | xmean1 = 1.-x1fac2/a ; |
---|
| 659 | |
---|
| 660 | const G4double fctail = factail*1.0 ; |
---|
| 661 | c = 2.+fctail*tau ; |
---|
| 662 | G4double c1 = c-1. ; |
---|
| 663 | G4double c2 = c-2. ; |
---|
| 664 | if(c2 == 0.) c2 = fctail*tausmall ; |
---|
| 665 | |
---|
| 666 | b = 1.+(c-xsi)/a ; |
---|
| 667 | |
---|
| 668 | b1 = b+1. ; |
---|
| 669 | bx = c/a ; |
---|
| 670 | eb1=exp((c1)*log(b1)) ; |
---|
| 671 | ebx=exp((c1)*log(bx)) ; |
---|
| 672 | xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/c2)/(eb1-ebx) ; |
---|
| 673 | |
---|
| 674 | f1x0 = a*ea/eaa ; |
---|
| 675 | f2x0 = c1*eb1*ebx/(eb1-ebx)/ |
---|
| 676 | exp(c*log(bx)) ; |
---|
| 677 | // from continuity at x=x0 |
---|
| 678 | prob = f2x0/(f1x0+f2x0) ; |
---|
| 679 | // from xmean = xmeanth |
---|
| 680 | qprob = (f1x0+f2x0)*xmeanth/(f2x0*xmean1+f1x0*xmean2) ; |
---|
| 681 | |
---|
| 682 | } |
---|
| 683 | |
---|
| 684 | // sampling of costheta |
---|
| 685 | if (G4UniformRand() < qprob) |
---|
| 686 | { |
---|
| 687 | if (G4UniformRand() < prob) |
---|
| 688 | cth = 1.+log(ea+G4UniformRand()*eaa)/a ; |
---|
| 689 | else |
---|
| 690 | cth = b-b1*bx/exp(log(ebx-G4UniformRand()*(ebx-eb1))/(c-1.)) ; |
---|
| 691 | } |
---|
| 692 | else |
---|
| 693 | { |
---|
| 694 | cth = -1.+2.*G4UniformRand(); |
---|
| 695 | } |
---|
| 696 | } |
---|
| 697 | } |
---|
| 698 | |
---|
| 699 | return cth ; |
---|
| 700 | } |
---|
| 701 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 702 | |
---|
| 703 | G4double G4MscModel71::SampleDisplacement() |
---|
| 704 | { |
---|
| 705 | const G4double kappa = 2.5; |
---|
| 706 | const G4double kappapl1 = kappa+1.; |
---|
| 707 | const G4double kappami1 = kappa-1.; |
---|
| 708 | G4double rmean = 0.0; |
---|
| 709 | if (currentTau >= tausmall) { |
---|
| 710 | if (currentTau < taulim) { |
---|
| 711 | rmean = kappa*currentTau*currentTau*currentTau*(1.-kappapl1*currentTau*0.25)/6. ; |
---|
| 712 | |
---|
| 713 | } else { |
---|
| 714 | G4double etau = 0.0; |
---|
| 715 | if (currentTau<taubig) etau = exp(-currentTau); |
---|
| 716 | rmean = -kappa*currentTau; |
---|
| 717 | rmean = -exp(rmean)/(kappa*kappami1); |
---|
| 718 | rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1; |
---|
| 719 | } |
---|
| 720 | if (rmean>0.) rmean = 2.*lambda0*sqrt(rmean/3.0); |
---|
| 721 | else rmean = 0.; |
---|
| 722 | } |
---|
| 723 | return rmean; |
---|
| 724 | } |
---|
| 725 | |
---|
| 726 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 727 | |
---|
| 728 | |
---|
| 729 | |
---|
| 730 | |
---|