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