CMS 3D CMS Logo

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