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