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