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

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

update ti head

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