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