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