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