source: trunk/source/processes/electromagnetic/standard/src/G4UrbanMscModel92.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

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