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