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