Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UrbanMscModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30// File name: G4UrbanMscModel
31//
32// Author: Laszlo Urban
33//
34// Creation date: 19.02.2013
35//
36// Created from G4UrbanMscModel96
37//
38// New parametrization for theta0
39// Correction for very small step length
40//
41// Class Description:
42//
43// Implementation of the model of multiple scattering based on
44// H.W.Lewis Phys Rev 78 (1950) 526 and others
45
46// -------------------------------------------------------------------
47// In its present form the model can be used for simulation
48// of the e-/e+ multiple scattering
49
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54#include "G4UrbanMscModel.hh"
56#include "G4SystemOfUnits.hh"
57#include "Randomize.hh"
58#include "G4Positron.hh"
59#include "G4EmParameters.hh"
62
63#include "G4Poisson.hh"
64#include "G4Pow.hh"
65#include "G4Log.hh"
66#include "G4Exp.hh"
67#include "G4AutoLock.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71std::vector<G4UrbanMscModel::mscData*> G4UrbanMscModel::msc;
72
73namespace
74{
75 G4Mutex theUrbanMutex = G4MUTEX_INITIALIZER;
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
81 : G4VMscModel(nam)
82{
83 masslimite = 0.6*CLHEP::MeV;
84 fr = 0.02;
85 taubig = 8.0;
86 tausmall = 1.e-16;
87 taulim = 1.e-6;
88 currentTau = taulim;
89 tlimitminfix = 0.01*CLHEP::nm;
90 tlimitminfix2 = 1.*CLHEP::nm;
91 stepmin = tlimitminfix;
92 smallstep = 1.e10;
93 currentRange = 0.;
94 rangeinit = 0.;
95 tlimit = 1.e10*CLHEP::mm;
96 tlimitmin = 10.*tlimitminfix;
97 tgeom = 1.e50*CLHEP::mm;
98 geombig = tgeom;
99 geommin = 1.e-6*CLHEP::mm;
100 geomlimit = geombig;
101 presafety = 0.;
102
103 positron = G4Positron::Positron();
104 rndmEngineMod = G4Random::getTheEngine();
105
106 drr = 0.35;
107 finalr = 10.*CLHEP::um;
108
109 tlow = 5.*CLHEP::keV;
110 invmev = 1.0/CLHEP::MeV;
111
112 skindepth = skin*stepmin;
113
114 mass = CLHEP::proton_mass_c2;
115 charge = chargeSquare = 1.0;
116 currentKinEnergy = currentRadLength = lambda0 = lambdaeff = tPathLength
117 = zPathLength = par1 = par2 = par3 = rndmarray[0] = rndmarray[1] = 0;
118 currentLogKinEnergy = LOG_EKIN_MIN;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
124{
125 if(isFirstInstance) {
126 for(auto const & ptr : msc) { delete ptr; }
127 msc.clear();
128 }
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
134 const G4DataVector&)
135{
136 // set values of some data members
137 SetParticle(p);
138 fParticleChange = GetParticleChangeForMSC(p);
140
141 latDisplasmentbackup = latDisplasment;
142
143 // if model is locked parameters should be defined via Set methods
144 if(!IsLocked()) {
146 fPosiCorrection = G4EmParameters::Instance()->MscPositronCorrection();
147 }
148
149 // initialise cache only once
150 if(0 == msc.size()) {
151 G4AutoLock l(&theUrbanMutex);
152 if(0 == msc.size()) {
153 isFirstInstance = true;
154 msc.resize(1, nullptr);
155 }
156 l.unlock();
157 }
158 // initialise cache for each new run
159 if(isFirstInstance) { InitialiseModelCache(); }
160
161 /*
162 G4cout << "### G4UrbanMscModel::Initialise done for "
163 << p->GetParticleName() << " type= " << steppingAlgorithm << G4endl;
164 G4cout << " RangeFact= " << facrange << " GeomFact= " << facgeom
165 << " SafetyFact= " << facsafety << " LambdaLim= " << lambdalimit
166 << G4endl;
167 */
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
173 const G4ParticleDefinition* part,
174 G4double kinEnergy,
175 G4double atomicNumber,G4double,
177{
178 static const G4double epsmin = 1.e-4 , epsmax = 1.e10;
179
180 static const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38.,47.,
181 50., 56., 64., 74., 79., 82. };
182
183 // corr. factors for e-/e+ lambda for T <= Tlim
184 static const G4double celectron[15][22] =
185 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
186 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
187 1.112,1.108,1.100,1.093,1.089,1.087 },
188 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
189 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
190 1.109,1.105,1.097,1.090,1.086,1.082 },
191 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
192 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
193 1.131,1.124,1.113,1.104,1.099,1.098 },
194 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
195 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
196 1.112,1.105,1.096,1.089,1.085,1.098 },
197 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
198 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
199 1.073,1.070,1.064,1.059,1.056,1.056 },
200 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
201 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
202 1.074,1.070,1.063,1.059,1.056,1.052 },
203 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
204 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
205 1.068,1.064,1.059,1.054,1.051,1.050 },
206 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
207 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
208 1.039,1.037,1.034,1.031,1.030,1.036 },
209 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
210 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
211 1.031,1.028,1.024,1.022,1.021,1.024 },
212 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
213 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
214 1.020,1.017,1.015,1.013,1.013,1.020 },
215 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
216 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
217 0.995,0.993,0.993,0.993,0.993,1.011 },
218 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
219 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
220 0.974,0.972,0.973,0.974,0.975,0.987 },
221 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
222 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
223 0.950,0.947,0.949,0.952,0.954,0.963 },
224 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
225 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
226 0.941,0.938,0.940,0.944,0.946,0.954 },
227 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
228 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
229 0.933,0.930,0.933,0.936,0.939,0.949 }};
230
231 static const G4double cpositron[15][22] = {
232 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
233 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
234 1.131,1.126,1.117,1.108,1.103,1.100 },
235 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
236 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
237 1.138,1.132,1.122,1.113,1.108,1.102 },
238 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
239 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
240 1.203,1.190,1.173,1.159,1.151,1.145 },
241 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
242 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
243 1.225,1.210,1.191,1.175,1.166,1.174 },
244 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
245 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
246 1.217,1.203,1.184,1.169,1.160,1.151 },
247 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
248 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
249 1.237,1.222,1.201,1.184,1.174,1.159 },
250 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
251 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
252 1.252,1.234,1.212,1.194,1.183,1.170 },
253 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
254 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
255 1.254,1.237,1.214,1.195,1.185,1.179 },
256 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
257 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
258 1.312,1.288,1.258,1.235,1.221,1.205 },
259 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
260 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
261 1.320,1.294,1.264,1.240,1.226,1.214 },
262 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
263 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
264 1.328,1.302,1.270,1.245,1.231,1.233 },
265 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
266 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
267 1.361,1.330,1.294,1.267,1.251,1.239 },
268 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
269 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
270 1.409,1.372,1.330,1.298,1.280,1.258 },
271 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
272 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
273 1.442,1.400,1.354,1.319,1.299,1.272 },
274 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
275 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
276 1.456,1.412,1.364,1.328,1.307,1.282 }};
277
278 //data/corrections for T > Tlim
279
280 static const G4double hecorr[15] = {
281 120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
282 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
283 -22.30};
284
285 G4double sigma;
286 SetParticle(part);
287
288 G4double Z23 = G4Pow::GetInstance()->Z23(G4lrint(atomicNumber));
289
290 // correction if particle .ne. e-/e+
291 // compute equivalent kinetic energy
292 // lambda depends on p*beta ....
293
294 G4double eKineticEnergy = kinEnergy;
295
296 if(mass > CLHEP::electron_mass_c2)
297 {
298 G4double TAU = kinEnergy/mass ;
299 G4double c = mass*TAU*(TAU+2.)/(CLHEP::electron_mass_c2*(TAU+1.)) ;
300 G4double w = c-2.;
301 G4double tau = 0.5*(w+std::sqrt(w*w+4.*c)) ;
302 eKineticEnergy = CLHEP::electron_mass_c2*tau ;
303 }
304
305 G4double eTotalEnergy = eKineticEnergy + CLHEP::electron_mass_c2 ;
306 G4double beta2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
307 /(eTotalEnergy*eTotalEnergy);
308 G4double bg2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
309 /(CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
310
311 static const G4double epsfactor = 2.*CLHEP::electron_mass_c2*
312 CLHEP::electron_mass_c2*CLHEP::Bohr_radius*CLHEP::Bohr_radius
313 /(CLHEP::hbarc*CLHEP::hbarc);
314 G4double eps = epsfactor*bg2/Z23;
315
316 if (eps<epsmin) sigma = 2.*eps*eps;
317 else if(eps<epsmax) sigma = G4Log(1.+2.*eps)-2.*eps/(1.+2.*eps);
318 else sigma = G4Log(2.*eps)-1.+1./eps;
319
320 sigma *= chargeSquare*atomicNumber*atomicNumber/(beta2*bg2);
321
322 // interpolate in AtomicNumber and beta2
323 G4double c1,c2,cc1;
324
325 // get bin number in Z
326 G4int iZ = 14;
327 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
328 while ((iZ>=0)&&(Zdat[iZ]>=atomicNumber)) { --iZ; }
329
330 iZ = std::min(std::max(iZ, 0), 13);
331
332 G4double ZZ1 = Zdat[iZ];
333 G4double ZZ2 = Zdat[iZ+1];
334 G4double ratZ = (atomicNumber-ZZ1)*(atomicNumber+ZZ1)/
335 ((ZZ2-ZZ1)*(ZZ2+ZZ1));
336
337 static const G4double Tlim = 10.*CLHEP::MeV;
338 static const G4double sigmafactor =
339 CLHEP::twopi*CLHEP::classic_electr_radius*CLHEP::classic_electr_radius;
340 static const G4double beta2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
341 ((Tlim+CLHEP::electron_mass_c2)*(Tlim+CLHEP::electron_mass_c2));
342 static const G4double bg2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
343 (CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
344
345 static const G4double sig0[15] = {
346 0.2672*CLHEP::barn, 0.5922*CLHEP::barn, 2.653*CLHEP::barn, 6.235*CLHEP::barn,
347 11.69*CLHEP::barn , 13.24*CLHEP::barn , 16.12*CLHEP::barn, 23.00*CLHEP::barn,
348 35.13*CLHEP::barn , 39.95*CLHEP::barn , 50.85*CLHEP::barn, 67.19*CLHEP::barn,
349 91.15*CLHEP::barn , 104.4*CLHEP::barn , 113.1*CLHEP::barn};
350
351 static const G4double Tdat[22] = {
352 100*CLHEP::eV, 200*CLHEP::eV, 400*CLHEP::eV, 700*CLHEP::eV,
353 1*CLHEP::keV, 2*CLHEP::keV, 4*CLHEP::keV, 7*CLHEP::keV,
354 10*CLHEP::keV, 20*CLHEP::keV, 40*CLHEP::keV, 70*CLHEP::keV,
355 100*CLHEP::keV, 200*CLHEP::keV, 400*CLHEP::keV, 700*CLHEP::keV,
356 1*CLHEP::MeV, 2*CLHEP::MeV, 4*CLHEP::MeV, 7*CLHEP::MeV,
357 10*CLHEP::MeV, 20*CLHEP::MeV};
358
359 if(eKineticEnergy <= Tlim)
360 {
361 // get bin number in T (beta2)
362 G4int iT = 21;
363 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
364 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
365
366 iT = std::min(std::max(iT, 0), 20);
367
368 // calculate betasquare values
369 G4double T = Tdat[iT];
370 G4double E = T + CLHEP::electron_mass_c2;
371 G4double b2small = T*(E+CLHEP::electron_mass_c2)/(E*E);
372
373 T = Tdat[iT+1];
374 E = T + CLHEP::electron_mass_c2;
375 G4double b2big = T*(E+CLHEP::electron_mass_c2)/(E*E);
376 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
377
378 if (charge < 0.)
379 {
380 c1 = celectron[iZ][iT];
381 c2 = celectron[iZ+1][iT];
382 cc1 = c1+ratZ*(c2-c1);
383
384 c1 = celectron[iZ][iT+1];
385 c2 = celectron[iZ+1][iT+1];
386 }
387 else
388 {
389 c1 = cpositron[iZ][iT];
390 c2 = cpositron[iZ+1][iT];
391 cc1 = c1+ratZ*(c2-c1);
392
393 c1 = cpositron[iZ][iT+1];
394 c2 = cpositron[iZ+1][iT+1];
395 }
396 G4double cc2 = c1+ratZ*(c2-c1);
397 sigma *= sigmafactor/(cc1+ratb2*(cc2-cc1));
398 }
399 else
400 {
401 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
402 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
403 if((atomicNumber >= ZZ1) && (atomicNumber <= ZZ2))
404 sigma = c1+ratZ*(c2-c1) ;
405 else if(atomicNumber < ZZ1)
406 sigma = atomicNumber*atomicNumber*c1/(ZZ1*ZZ1);
407 else if(atomicNumber > ZZ2)
408 sigma = atomicNumber*atomicNumber*c2/(ZZ2*ZZ2);
409 }
410 // low energy correction based on theory
411 sigma *= (1.+0.30/(1.+std::sqrt(1000.*eKineticEnergy)));
412
413 return sigma;
414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
417
419{
420 SetParticle(track->GetDynamicParticle()->GetDefinition());
421 firstStep = true;
422 insideskin = false;
423 fr = facrange;
424 tlimit = tgeom = rangeinit = geombig;
425 smallstep = 1.e10;
426 stepmin = tlimitminfix;
427 tlimitmin = 10.*tlimitminfix;
428 rndmEngineMod = G4Random::getTheEngine();
429}
430
431//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
432
434 const G4Track& track,
435 G4double& currentMinimalStep)
436{
437 tPathLength = currentMinimalStep;
438 const G4DynamicParticle* dp = track.GetDynamicParticle();
439
440 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
441 G4StepStatus stepStatus = sp->GetStepStatus();
442 couple = track.GetMaterialCutsCouple();
443 SetCurrentCouple(couple);
444 idx = couple->GetIndex();
445 currentKinEnergy = dp->GetKineticEnergy();
446 currentLogKinEnergy = dp->GetLogKineticEnergy();
447 currentRange = GetRange(particle,currentKinEnergy,couple,currentLogKinEnergy);
448 lambda0 = GetTransportMeanFreePath(particle,currentKinEnergy,
449 currentLogKinEnergy);
450 tPathLength = std::min(tPathLength,currentRange);
451 /*
452 G4cout << "G4Urban::StepLimit tPathLength= " << tPathLength
453 << " range= " <<currentRange<< " lambda= "<<lambda0
454 <<G4endl;
455 */
456 // extreme small step
457 if(tPathLength < tlimitminfix) {
458 latDisplasment = false;
459 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
460 }
461
462 presafety = (stepStatus == fGeomBoundary) ? sp->GetSafety()
463 : ComputeSafety(sp->GetPosition(), tPathLength);
464
465 // stop here if small step or range is less than safety
466 if(tPathLength == currentRange && tPathLength < presafety) {
467 latDisplasment = false;
468 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
469 }
470
471 // upper limit for the straight line distance the particle can travel
472 // for electrons and positrons
473 G4double distance = (mass < masslimite)
474 ? currentRange*msc[idx]->doverra
475 // for muons, hadrons
476 : currentRange*msc[idx]->doverrb;
477
478 /*
479 G4cout << "G4Urban::StepLimit tPathLength= "
480 <<tPathLength<<" safety= " << presafety
481 << " range= " <<currentRange<< " lambda= "<<lambda0
482 << " Alg: " << steppingAlgorithm <<G4endl;
483 */
484 // far from geometry boundary
485 if(distance < presafety)
486 {
487 latDisplasment = false;
488 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
489 }
490
491 latDisplasment = latDisplasmentbackup;
492 // ----------------------------------------------------------------
493 // distance to boundary
495 {
496 //compute geomlimit and presafety
497 geomlimit = ComputeGeomLimit(track, presafety, currentRange);
498 /*
499 G4cout << "G4Urban::Distance to boundary geomlimit= "
500 <<geomlimit<<" safety= " << presafety<<G4endl;
501 */
502
503 smallstep += 1.;
504 insideskin = false;
505 tgeom = geombig;
506
507 // initialisation at first step and at the boundary
508 if(firstStep || (stepStatus == fGeomBoundary))
509 {
510 rangeinit = currentRange;
511 if(!firstStep) { smallstep = 1.; }
512
513 //stepmin ~ lambda_elastic
514 stepmin = ComputeStepmin();
515 skindepth = skin*stepmin;
516 tlimitmin = ComputeTlimitmin();
517 /*
518 G4cout << "rangeinit= " << rangeinit << " stepmin= " << stepmin
519 << " tlimitmin= " << tlimitmin << " geomlimit= "
520 << geomlimit <<G4endl;
521 */
522 }
523
524 // constraint from the geometry
525 // tgeom is upper limit for the step size
526 if((geomlimit < geombig) && (geomlimit > geommin))
527 {
528 // geomlimit is a geometrical step length
529 // transform it to true path length (estimation)
530 if(lambda0 > geomlimit) {
531 geomlimit = -lambda0*G4Log(1.-geomlimit/lambda0)+tlimitmin;
532 }
533 tgeom = (stepStatus == fGeomBoundary) ? geomlimit/facgeom
534 : facrange*rangeinit + stepmin;
535 }
536 else if(geomlimit > geombig) {
537 // range smaller than distance to boundary
538 // here tgeom is the true path length
539 tgeom = currentRange;
540 }
541 else if((geomlimit < geommin) && (geomlimit > 0.)) {
542 // geomlimit small (smaller than geommin=1 um)
543 // here true pathlength ~ geom path length
544 tgeom = geomlimit;
545 }
546
547 //step limit
548 tlimit = (currentRange > presafety) ? facrange*rangeinit : currentRange;
549
550 //lower limit for tlimit
551 tlimit = std::max(tlimit, tlimitmin);
552 tlimit = std::min(tlimit, tgeom);
553 /*
554 G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit
555 << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
556 */
557 // shortcut
558 if((tPathLength < tlimit) && (tPathLength < presafety) &&
559 (smallstep > skin) && (tPathLength < geomlimit-0.999*skindepth))
560 {
561 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
562 }
563
564 // step reduction near to boundary
565 if(smallstep <= skin)
566 {
567 tlimit = stepmin;
568 insideskin = true;
569 }
570 else if(geomlimit < geombig)
571 {
572 if(geomlimit > skindepth)
573 {
574 tlimit = std::min(tlimit, geomlimit-0.999*skindepth);
575 }
576 else
577 {
578 insideskin = true;
579 tlimit = std::min(tlimit, stepmin);
580 }
581 }
582
583 tlimit = std::max(tlimit, stepmin);
584
585 // randomise if not 'small' step and step determined by msc
586 tPathLength = (tlimit < tPathLength && smallstep > skin && !insideskin)
587 ? std::min(tPathLength, Randomizetlimit())
588 : std::min(tPathLength, tlimit);
589 }
590 // ----------------------------------------------------------------
591 // for simulation with or without magnetic field
592 // there no small step/single scattering at boundaries
593 else if(steppingAlgorithm == fUseSafety)
594 {
595 if(firstStep || (stepStatus == fGeomBoundary)) {
596 rangeinit = currentRange;
597 fr = facrange;
598 // stepping for e+/e- only (not for muons,hadrons)
599 if(mass < masslimite)
600 {
601 rangeinit = std::max(rangeinit, lambda0);
602 if(lambda0 > lambdalimit) {
603 fr *= (0.75+0.25*lambda0/lambdalimit);
604 }
605 }
606 //lower limit for tlimit
607 stepmin = ComputeStepmin();
608 tlimitmin = ComputeTlimitmin();
609 }
610
611 //step limit
612 tlimit = (currentRange > presafety) ?
613 std::max(fr*rangeinit, facsafety*presafety) : currentRange;
614
615 //lower limit for tlimit
616 tlimit = std::max(tlimit, tlimitmin);
617
618 // randomise if step determined by msc
619 tPathLength = (tlimit < tPathLength) ?
620 std::min(tPathLength, Randomizetlimit()) : tPathLength;
621 }
622 // ----------------------------------------------------------------
623 // for simulation with or without magnetic field
624 // there is small step/single scattering at boundaries
626 {
627 if(firstStep || (stepStatus == fGeomBoundary)) {
628 rangeinit = currentRange;
629 fr = facrange;
630 if(mass < masslimite)
631 {
632 if(lambda0 > lambdalimit) {
633 fr *= (0.84+0.16*lambda0/lambdalimit);
634 }
635 }
636 //lower limit for tlimit
637 stepmin = ComputeStepmin();
638 tlimitmin = ComputeTlimitmin();
639 }
640 //step limit
641 tlimit = (currentRange > presafety) ?
642 std::max(fr*rangeinit, facsafety*presafety) : currentRange;
643
644 //lower limit for tlimit
645 tlimit = std::max(tlimit, tlimitmin);
646
647 // condition for tPathLength from drr and finalr
648 if(currentRange > finalr) {
649 G4double tmax = drr*currentRange+
650 finalr*(1.-drr)*(2.-finalr/currentRange);
651 tPathLength = std::min(tPathLength,tmax);
652 }
653
654 // randomise if step determined by msc
655 tPathLength = (tlimit < tPathLength) ?
656 std::min(tPathLength, Randomizetlimit()) : tPathLength;
657 }
658
659 // ----------------------------------------------------------------
660 // simple step limitation
661 else
662 {
663 if (stepStatus == fGeomBoundary)
664 {
665 tlimit = (currentRange > lambda0)
666 ? facrange*currentRange : facrange*lambda0;
667 tlimit = std::max(tlimit, tlimitmin);
668 }
669 // randomise if step determined by msc
670 tPathLength = (tlimit < tPathLength) ?
671 std::min(tPathLength, Randomizetlimit()) : tPathLength;
672 }
673
674 // ----------------------------------------------------------------
675 firstStep = false;
676 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
677}
678
679//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
680
682{
683 lambdaeff = lambda0;
684 par1 = -1. ;
685 par2 = par3 = 0. ;
686
687 // this correction needed to run MSC with eIoni and eBrem inactivated
688 // and makes no harm for a normal run
689 tPathLength = std::min(tPathLength,currentRange);
690
691 // do the true -> geom transformation
692 zPathLength = tPathLength;
693
694 // z = t for very small tPathLength
695 if(tPathLength < tlimitminfix2) return zPathLength;
696
697 /*
698 G4cout << "ComputeGeomPathLength: tpl= " << tPathLength
699 << " R= " << currentRange << " L0= " << lambda0
700 << " E= " << currentKinEnergy << " "
701 << particle->GetParticleName() << G4endl;
702 */
703 G4double tau = tPathLength/lambda0 ;
704
705 if ((tau <= tausmall) || insideskin) {
706 zPathLength = std::min(tPathLength, lambda0);
707
708 } else if (tPathLength < currentRange*dtrl) {
709 zPathLength = (tau < taulim) ? tPathLength*(1.-0.5*tau)
710 : lambda0*(1.-G4Exp(-tau));
711
712 } else if(currentKinEnergy < mass || tPathLength == currentRange) {
713 par1 = 1./currentRange;
714 par2 = currentRange/lambda0;
715 par3 = 1.+par2;
716 if(tPathLength < currentRange) {
717 zPathLength =
718 (1.-G4Exp(par3*G4Log(1.-tPathLength/currentRange)))/(par1*par3);
719 } else {
720 zPathLength = 1./(par1*par3);
721 }
722
723 } else {
724 G4double rfin = std::max(currentRange-tPathLength, 0.01*currentRange);
725 G4double T1 = GetEnergy(particle,rfin,couple);
726 G4double lambda1 = GetTransportMeanFreePath(particle,T1);
727
728 par1 = (lambda0-lambda1)/(lambda0*tPathLength);
729 //G4cout << "par1= " << par1 << " L1= " << lambda1 << G4endl;
730 par2 = 1./(par1*lambda0);
731 par3 = 1.+par2;
732 zPathLength = (1.-G4Exp(par3*G4Log(lambda1/lambda0)))/(par1*par3);
733 }
734
735 zPathLength = std::min(zPathLength, lambda0);
736 //G4cout<< "zPathLength= "<< zPathLength<< " L0= " << lambda0 << G4endl;
737 return zPathLength;
738}
739
740//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
741
743{
744 // step defined other than transportation
745 if(geomStepLength == zPathLength) {
746 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
747 // << " step= " << geomStepLength << " *** " << G4endl;
748 return tPathLength;
749 }
750
751 zPathLength = geomStepLength;
752
753 // t = z for very small step
754 if(geomStepLength < tlimitminfix2) {
755 tPathLength = geomStepLength;
756
757 // recalculation
758 } else {
759
760 G4double tlength = geomStepLength;
761 if((geomStepLength > lambda0*tausmall) && !insideskin) {
762
763 if(par1 < 0.) {
764 tlength = -lambda0*G4Log(1.-geomStepLength/lambda0) ;
765 } else {
766 const G4double par4 = par1*par3;
767 if(par4*geomStepLength < 1.) {
768 tlength = (1.-G4Exp(G4Log(1.-par4*geomStepLength)/par3))/par1;
769 } else {
770 tlength = currentRange;
771 }
772 }
773
774 if(tlength < geomStepLength) { tlength = geomStepLength; }
775 else if(tlength > tPathLength) { tlength = tPathLength; }
776 }
777 tPathLength = tlength;
778 }
779 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
780 // << " step= " << geomStepLength << " &&& " << G4endl;
781
782 return tPathLength;
783}
784
785//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
786
789 G4double /*safety*/)
790{
791 fDisplacement.set(0.0,0.0,0.0);
792 if(tPathLength >= currentRange) { return fDisplacement; }
793
794 G4double kinEnergy = currentKinEnergy;
795 if (tPathLength > currentRange*dtrl) {
796 kinEnergy = GetEnergy(particle,currentRange-tPathLength,couple);
797 } else if(tPathLength > currentRange*0.01) {
798 kinEnergy -= tPathLength*GetDEDX(particle,currentKinEnergy,couple,
799 currentLogKinEnergy);
800 }
801
802 if((tPathLength <= tlimitminfix) || (tPathLength < tausmall*lambda0) ||
803 (kinEnergy <= CLHEP::eV)) { return fDisplacement; }
804
805 G4double cth = SampleCosineTheta(tPathLength,kinEnergy);
806
807 // protection against 'bad' cth values
808 if(std::abs(cth) >= 1.0) { return fDisplacement; }
809
810 G4double sth = std::sqrt((1.0 - cth)*(1.0 + cth));
811 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
812 G4ThreeVector newDirection(sth*std::cos(phi),sth*std::sin(phi),cth);
813 newDirection.rotateUz(oldDirection);
814
815 fParticleChange->ProposeMomentumDirection(newDirection);
816 /*
817 G4cout << "G4UrbanMscModel::SampleSecondaries: e(MeV)= " << kineticEnergy
818 << " sinTheta= " << sth << " safety(mm)= " << safety
819 << " trueStep(mm)= " << tPathLength
820 << " geomStep(mm)= " << zPathLength
821 << G4endl;
822 */
823
824 if (latDisplasment && currentTau >= tausmall) {
825 if(dispAlg96) { SampleDisplacement(sth, phi); }
826 else { SampleDisplacementNew(cth, phi); }
827 fDisplacement.rotateUz(oldDirection);
828 }
829 return fDisplacement;
830}
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
833
834G4double G4UrbanMscModel::SampleCosineTheta(G4double trueStepLength,
835 G4double kinEnergy)
836{
837 G4double cth = 1.0;
838 G4double tau = trueStepLength/lambda0;
839
840 // mean tau value
841 if(currentKinEnergy != kinEnergy) {
842 G4double lambda1 = GetTransportMeanFreePath(particle, kinEnergy);
843 if(std::abs(lambda1 - lambda0) > lambda0*0.01 && lambda1 > 0.) {
844 tau = trueStepLength*G4Log(lambda0/lambda1)/(lambda0-lambda1);
845 }
846 }
847
848 currentTau = tau;
849 lambdaeff = trueStepLength/currentTau;
850 currentRadLength = couple->GetMaterial()->GetRadlen();
851
852 if (tau >= taubig) { cth = -1.+2.*rndmEngineMod->flat(); }
853 else if (tau >= tausmall) {
854 static const G4double numlim = 0.01;
855 static const G4double onethird = 1./3.;
856 if(tau < numlim) {
857 xmeanth = 1.0 - tau*(1.0 - 0.5*tau);
858 x2meanth= 1.0 - tau*(5.0 - 6.25*tau)*onethird;
859 } else {
860 xmeanth = G4Exp(-tau);
861 x2meanth = (1.+2.*G4Exp(-2.5*tau))*onethird;
862 }
863
864 // too large step of low-energy particle
865 G4double relloss = 1. - kinEnergy/currentKinEnergy;
866 static const G4double rellossmax= 0.50;
867 if(relloss > rellossmax) {
868 return SimpleScattering();
869 }
870 // is step extreme small ?
871 G4bool extremesmallstep = false;
872 G4double tsmall = std::min(tlimitmin,lambdalimit);
873
874 G4double theta0;
875 if(trueStepLength > tsmall) {
876 theta0 = ComputeTheta0(trueStepLength,kinEnergy);
877 } else {
878 theta0 = std::sqrt(trueStepLength/tsmall)
879 *ComputeTheta0(tsmall,kinEnergy);
880 extremesmallstep = true;
881 }
882
883 static const G4double onesixth = 1./6.;
884 static const G4double one12th = 1./12.;
885 static const G4double theta0max = CLHEP::pi*onesixth;
886
887 // protection for very small angles
888 G4double theta2 = theta0*theta0;
889
890 if(theta2 < tausmall) { return cth; }
891 if(theta0 > theta0max) { return SimpleScattering(); }
892
893 G4double x = theta2*(1.0 - theta2*one12th);
894 if(theta2 > numlim) {
895 G4double sth = 2*std::sin(0.5*theta0);
896 x = sth*sth;
897 }
898
899 // parameter for tail
900 G4double ltau = G4Log(tau);
901 G4double u = !extremesmallstep ? G4Exp(ltau*onesixth)
902 : G4Exp(G4Log(tsmall/lambda0)*onesixth);
903
904 G4double xx = G4Log(lambdaeff/currentRadLength);
905 G4double xsi = msc[idx]->coeffc1 +
906 u*(msc[idx]->coeffc2+msc[idx]->coeffc3*u)+msc[idx]->coeffc4*xx;
907
908 // tail should not be too big
909 xsi = std::max(xsi, 1.9);
910 /*
911 if(KineticEnergy > 20*MeV && xsi < 1.6) {
912 G4cout << "G4UrbanMscModel::SampleCosineTheta: E(GeV)= "
913 << KineticEnergy/GeV
914 << " !!** c= " << xsi
915 << " **!! length(mm)= " << trueStepLength << " Zeff= " << Zeff
916 << " " << couple->GetMaterial()->GetName()
917 << " tau= " << tau << G4endl;
918 }
919 */
920
921 G4double c = xsi;
922
923 if(std::abs(c-3.) < 0.001) { c = 3.001; }
924 else if(std::abs(c-2.) < 0.001) { c = 2.001; }
925
926 G4double c1 = c-1.;
927 G4double ea = G4Exp(-xsi);
928 G4double eaa = 1.-ea ;
929 G4double xmean1 = 1.-(1.-(1.+xsi)*ea)*x/eaa;
930 G4double x0 = 1. - xsi*x;
931
932 // G4cout << " xmean1= " << xmean1 << " xmeanth= " << xmeanth << G4endl;
933
934 if(xmean1 <= 0.999*xmeanth) { return SimpleScattering(); }
935
936 //from continuity of derivatives
937 G4double b = 1.+(c-xsi)*x;
938
939 G4double b1 = b+1.;
940 G4double bx = c*x;
941
942 G4double eb1 = G4Exp(G4Log(b1)*c1);
943 G4double ebx = G4Exp(G4Log(bx)*c1);
944 G4double d = ebx/eb1;
945
946 G4double xmean2 = (x0 + d - (bx - b1*d)/(c-2.))/(1. - d);
947
948 G4double f1x0 = ea/eaa;
949 G4double f2x0 = c1/(c*(1. - d));
950 G4double prob = f2x0/(f1x0+f2x0);
951
952 G4double qprob = xmeanth/(prob*xmean1+(1.-prob)*xmean2);
953
954 // sampling of costheta
955 //G4cout << "c= " << c << " qprob= " << qprob << " eb1= " << eb1
956 // << " c1= " << c1 << " b1= " << b1 << " bx= " << bx << " eb1= " << eb1
957 // << G4endl;
958 rndmEngineMod->flatArray(2, rndmarray);
959 if(rndmarray[0] < qprob)
960 {
961 G4double var = 0;
962 if(rndmarray[1] < prob) {
963 cth = 1.+G4Log(ea+rndmEngineMod->flat()*eaa)*x;
964 } else {
965 var = (1.0 - d)*rndmEngineMod->flat();
966 if(var < numlim*d) {
967 var /= (d*c1);
968 cth = -1.0 + var*(1.0 - 0.5*var*c)*(2. + (c - xsi)*x);
969 } else {
970 cth = 1. + x*(c - xsi - c*G4Exp(-G4Log(var + d)/c1));
971 }
972 }
973 } else {
974 cth = -1.+2.*rndmarray[1];
975 }
976 }
977 return cth;
978}
979
980//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
981
983 G4double kinEnergy)
984{
985 // for all particles take the width of the central part
986 // from a parametrization similar to the Highland formula
987 // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
988 G4double invbetacp = (kinEnergy+mass)/(kinEnergy*(kinEnergy+2.*mass));
989 if(currentKinEnergy != kinEnergy) {
990 invbetacp = std::sqrt(invbetacp*(currentKinEnergy+mass)/
991 (currentKinEnergy*(currentKinEnergy+2.*mass)));
992 }
993 G4double y = trueStepLength/currentRadLength;
994
995 if(fPosiCorrection && particle == positron)
996 {
997 static const G4double xl= 0.6;
998 static const G4double xh= 0.9;
999 static const G4double e = 113.0;
1000 G4double corr;
1001
1002 G4double tau = std::sqrt(currentKinEnergy*kinEnergy)/mass;
1003 G4double x = std::sqrt(tau*(tau+2.)/((tau+1.)*(tau+1.)));
1004 G4double a = msc[idx]->posa;
1005 G4double b = msc[idx]->posb;
1006 G4double c = msc[idx]->posc;
1007 G4double d = msc[idx]->posd;
1008 if(x < xl) {
1009 corr = a*(1.-G4Exp(-b*x));
1010 } else if(x > xh) {
1011 corr = c+d*G4Exp(e*(x-1.));
1012 } else {
1013 G4double yl = a*(1.-G4Exp(-b*xl));
1014 G4double yh = c+d*G4Exp(e*(xh-1.));
1015 G4double y0 = (yh-yl)/(xh-xl);
1016 G4double y1 = yl-y0*xl;
1017 corr = y0*x+y1;
1018 }
1019 //==================================================================
1020 y *= corr*msc[idx]->pose;
1021 }
1022
1023 static const G4double c_highland = 13.6*CLHEP::MeV;
1024 G4double theta0 = c_highland*std::abs(charge)*std::sqrt(y)*invbetacp;
1025
1026 // correction factor from e- scattering data
1027 theta0 *= (msc[idx]->coeffth1+msc[idx]->coeffth2*G4Log(y));
1028 return theta0;
1029}
1030
1031//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1032
1033void G4UrbanMscModel::SampleDisplacement(G4double, G4double phi)
1034{
1035 // simple and fast sampling
1036 // based on single scattering results
1037 // u = r/rmax : mean value
1038
1039 G4double rmax = std::sqrt((tPathLength-zPathLength)*(tPathLength+zPathLength));
1040 if(rmax > 0.)
1041 {
1042 G4double r = 0.73*rmax;
1043
1044 // simple distribution for v=Phi-phi=psi ~exp(-beta*v)
1045 // beta determined from the requirement that distribution should give
1046 // the same mean value than that obtained from the ss simulation
1047
1048 static const G4double cbeta = 2.160;
1049 static const G4double cbeta1 = 1. - G4Exp(-cbeta*CLHEP::pi);
1050 rndmEngineMod->flatArray(2, rndmarray);
1051 G4double psi = -G4Log(1. - rndmarray[0]*cbeta1)/cbeta;
1052 G4double Phi = (rndmarray[1] < 0.5) ? phi+psi : phi-psi;
1053 fDisplacement.set(r*std::cos(Phi),r*std::sin(Phi),0.0);
1054 }
1055}
1056
1057//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1058
1059void G4UrbanMscModel::SampleDisplacementNew(G4double, G4double phi)
1060{
1061 // simple and fast sampling
1062 // based on single scattering results
1063 // u = (r/rmax)**2 : distribution from ss simulations
1064 const G4double eps = 1.e-3;
1065 const G4double rmax =
1066 std::sqrt((tPathLength-zPathLength)*(tPathLength+zPathLength));
1067
1068 if(rmax > 0.)
1069 {
1070 const G4double x0 = 0.73 ;
1071 const G4double alpha = G4Log(7.33)/x0 ;
1072 const G4double a1 = 1.-x0 ;
1073 const G4double a2 = 1.-G4Exp(-alpha*x0) ;
1074 const G4double a3 = G4Exp(alpha*x0)-1. ;
1075 const G4double w1 = 2.*a2/(alpha*a1+2.*a2) ;
1076
1077 G4double r, sqx;
1078 if (rmax/currentRange < eps)
1079 {
1080 r = 0.73*rmax ;
1081 sqx = 1.;
1082 }
1083 else
1084 {
1085 rndmEngineMod->flatArray(2,rndmarray);
1086 const G4double x = (rndmarray[0] < w1) ? G4Log(1. + a3*rndmarray[1])/alpha :
1087 1. - a1*std::sqrt(1.-rndmarray[1]);
1088
1089 sqx = std::sqrt(x);
1090 r = sqx*rmax;
1091 }
1092 // Gaussian distribution for Phi-phi=psi
1093 const G4double sigma = 0.1+0.9*sqx;
1094 const G4double psi = G4RandGauss::shoot(0.,sigma);
1095 const G4double Phi = phi+psi;
1096 fDisplacement.set(r*std::cos(Phi), r*std::sin(Phi), 0.0);
1097 }
1098}
1099
1100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1101
1102void G4UrbanMscModel::InitialiseModelCache()
1103{
1104 // it is assumed, that for the second run only addition
1105 // of a new G4MaterialCutsCouple is possible
1106 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
1107 std::size_t numOfCouples = theCoupleTable->GetTableSize();
1108 if(numOfCouples != msc.size()) { msc.resize(numOfCouples, nullptr); }
1109
1110 for(G4int j=0; j<(G4int)numOfCouples; ++j) {
1111 auto aCouple = theCoupleTable->GetMaterialCutsCouple(j);
1112
1113 // new couple
1114 msc[j] = new mscData();
1115 G4double Zeff = aCouple->GetMaterial()->GetIonisation()->GetZeffective();
1116 G4double sqrz = std::sqrt(Zeff);
1117 msc[j]->sqrtZ = sqrz;
1118 // parameterisation of step limitation
1119 msc[j]->factmin = dispAlg96 ? 0.001 : 0.001/(1.+0.028*sqrz);
1120 G4double lnZ = G4Log(Zeff);
1121 // correction in theta0 formula
1122 G4double w = G4Exp(lnZ/6.);
1123 G4double facz = 0.990395+w*(-0.168386+w*0.093286);
1124 msc[j]->coeffth1 = facz*(1. - 8.7780e-2/Zeff);
1125 msc[j]->coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
1126
1127 // tail parameters
1128 G4double Z13 = w*w;
1129 msc[j]->coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
1130 msc[j]->coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
1131 msc[j]->coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
1132 msc[j]->coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
1133
1134 msc[j]->Z23 = Z13*Z13;
1135
1136 msc[j]->stepmina = 27.725/(1.+0.203*Zeff);
1137 msc[j]->stepminb = 6.152/(1.+0.111*Zeff);
1138
1139 // 21.07.2020
1140 msc[j]->doverra = 9.6280e-1 - 8.4848e-2*msc[j]->sqrtZ + 4.3769e-3*Zeff;
1141
1142 // 06.10.2020
1143 // msc[j]->doverra = 7.7024e-1 - 6.7878e-2*msc[j]->sqrtZ + 3.5015e-3*Zeff;
1144 msc[j]->doverrb = 1.15 - 9.76e-4*Zeff;
1145
1146 // corrections for e+
1147 msc[j]->posa = 0.994-4.08e-3*Zeff;
1148 msc[j]->posb = 7.16+(52.6+365./Zeff)/Zeff;
1149 msc[j]->posc = 1.000-4.47e-3*Zeff;
1150 msc[j]->posd = 1.21e-3*Zeff;
1151 msc[j]->pose = 1.+Zeff*(1.84035e-4*Zeff-1.86427e-2)+0.41125;
1152 }
1153}
1154
1155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
G4StepStatus
@ fGeomBoundary
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
virtual void flatArray(const int size, double *vect)=0
G4double GetLogKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4bool LateralDisplacementAlg96() const
static G4EmParameters * Instance()
G4bool MscPositronCorrection() const
static G4Positron * Positron()
Definition G4Positron.cc:90
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double Z23(G4int Z) const
Definition G4Pow.hh:125
static G4ProductionCutsTable * GetProductionCutsTable()
G4StepPoint * GetPreStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double ComputeTrueStepLength(G4double geomStepLength) override
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
void StartTracking(G4Track *) override
~G4UrbanMscModel() override
G4UrbanMscModel(const G4String &nam="UrbanMsc")
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeGeomPathLength(G4double truePathLength) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4bool IsLocked() const
G4double dtrl
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double facrange
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
G4double skin
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
G4VMscModel(const G4String &nam)
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double lambdalimit
G4MscStepLimitType steppingAlgorithm
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
G4bool latDisplasment
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
G4double facsafety
G4ThreeVector fDisplacement
void InitialiseParameters(const G4ParticleDefinition *)
G4double facgeom
#define LOG_EKIN_MIN
Definition templates.hh:98
int G4lrint(double ad)
Definition templates.hh:134