CMS 3D CMS Logo

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