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