Geant4-11
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//
29// GEANT4 Class file
30//
31//
32// File name: G4UrbanMscModel
33//
34// Author: Laszlo Urban
35//
36// Creation date: 19.02.2013
37//
38// Created from G4UrbanMscModel96
39//
40// New parametrization for theta0
41// Correction for very small step length
42//
43// Class Description:
44//
45// Implementation of the model of multiple scattering based on
46// H.W.Lewis Phys Rev 78 (1950) 526 and others
47
48// -------------------------------------------------------------------
49// In its present form the model can be used for simulation
50// of the e-/e+ multiple scattering
51//
52
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
57#include "G4UrbanMscModel.hh"
59#include "G4SystemOfUnits.hh"
60#include "Randomize.hh"
61#include "G4Electron.hh"
62#include "G4Positron.hh"
63#include "G4LossTableManager.hh"
64#include "G4EmParameters.hh"
67
68#include "G4Poisson.hh"
69#include "G4Pow.hh"
70#include "globals.hh"
71#include "G4Log.hh"
72#include "G4Exp.hh"
73#include "G4Threading.hh"
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
77std::vector<G4UrbanMscModel::mscData*> G4UrbanMscModel::msc;
78
79namespace
80{
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87 : G4VMscModel(nam)
88{
90 fr = 0.02;
91 taubig = 8.0;
92 tausmall = 1.e-16;
93 taulim = 1.e-6;
95 tlimitminfix = 0.01*CLHEP::nm;
98 smallstep = 1.e10;
99 currentRange = 0. ;
100 rangeinit = 0.;
101 tlimit = 1.e10*CLHEP::mm;
102 tlimitmin = 10.*tlimitminfix;
103 tgeom = 1.e50*CLHEP::mm;
104 geombig = tgeom;
105 geommin = 1.e-3*CLHEP::mm;
107 presafety = 0.*CLHEP::mm;
108
109 particle = nullptr;
110
113 rndmEngineMod = G4Random::getTheEngine();
114
115 firstStep = true;
116 insideskin = false;
117 latDisplasmentbackup = false;
118 dispAlg96 = true;
119
121 drr = 0.35;
122 finalr = 10.*CLHEP::um;
123
124 tlow = 5.*CLHEP::keV;
125 invmev = 1.0/CLHEP::MeV;
126
128
130 charge = chargeSquare = 1.0;
132 = zPathLength = par1 = par2 = par3 = rndmarray[0] = rndmarray[1] = 0;
134
135 idx = 0;
136 fParticleChange = nullptr;
137 couple = nullptr;
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
143{
144 if(isFirstInstance) {
145 for(auto & ptr : msc) { delete ptr; }
146 msc.clear();
147 }
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
153 const G4DataVector&)
154{
155 // set values of some data members
156 SetParticle(p);
159
162
163 // initialise cache only once
164 if(0 == msc.size()) {
166 if(0 == msc.size()) {
167 isFirstInstance = true;
168 msc.resize(1, nullptr);
169 }
170 l.unlock();
171 }
172 // initialise cache for each new run
174
175 /*
176 G4cout << "### G4UrbanMscModel::Initialise done for "
177 << p->GetParticleName() << " type= " << steppingAlgorithm << G4endl;
178 G4cout << " RangeFact= " << facrange << " GeomFact= " << facgeom
179 << " SafetyFact= " << facsafety << " LambdaLim= " << lambdalimit
180 << G4endl;
181 */
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185
187 const G4ParticleDefinition* part,
188 G4double kinEnergy,
189 G4double atomicNumber,G4double,
191{
192 static const G4double epsmin = 1.e-4 , epsmax = 1.e10;
193
194 static const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38.,47.,
195 50., 56., 64., 74., 79., 82. };
196
197 // corr. factors for e-/e+ lambda for T <= Tlim
198 static const G4double celectron[15][22] =
199 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
200 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
201 1.112,1.108,1.100,1.093,1.089,1.087 },
202 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
203 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
204 1.109,1.105,1.097,1.090,1.086,1.082 },
205 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
206 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
207 1.131,1.124,1.113,1.104,1.099,1.098 },
208 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
209 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
210 1.112,1.105,1.096,1.089,1.085,1.098 },
211 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
212 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
213 1.073,1.070,1.064,1.059,1.056,1.056 },
214 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
215 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
216 1.074,1.070,1.063,1.059,1.056,1.052 },
217 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
218 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
219 1.068,1.064,1.059,1.054,1.051,1.050 },
220 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
221 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
222 1.039,1.037,1.034,1.031,1.030,1.036 },
223 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
224 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
225 1.031,1.028,1.024,1.022,1.021,1.024 },
226 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
227 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
228 1.020,1.017,1.015,1.013,1.013,1.020 },
229 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
230 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
231 0.995,0.993,0.993,0.993,0.993,1.011 },
232 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
233 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
234 0.974,0.972,0.973,0.974,0.975,0.987 },
235 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
236 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
237 0.950,0.947,0.949,0.952,0.954,0.963 },
238 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
239 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
240 0.941,0.938,0.940,0.944,0.946,0.954 },
241 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
242 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
243 0.933,0.930,0.933,0.936,0.939,0.949 }};
244
245 static const G4double cpositron[15][22] = {
246 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
247 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
248 1.131,1.126,1.117,1.108,1.103,1.100 },
249 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
250 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
251 1.138,1.132,1.122,1.113,1.108,1.102 },
252 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
253 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
254 1.203,1.190,1.173,1.159,1.151,1.145 },
255 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
256 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
257 1.225,1.210,1.191,1.175,1.166,1.174 },
258 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
259 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
260 1.217,1.203,1.184,1.169,1.160,1.151 },
261 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
262 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
263 1.237,1.222,1.201,1.184,1.174,1.159 },
264 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
265 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
266 1.252,1.234,1.212,1.194,1.183,1.170 },
267 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
268 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
269 1.254,1.237,1.214,1.195,1.185,1.179 },
270 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
271 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
272 1.312,1.288,1.258,1.235,1.221,1.205 },
273 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
274 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
275 1.320,1.294,1.264,1.240,1.226,1.214 },
276 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
277 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
278 1.328,1.302,1.270,1.245,1.231,1.233 },
279 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
280 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
281 1.361,1.330,1.294,1.267,1.251,1.239 },
282 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
283 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
284 1.409,1.372,1.330,1.298,1.280,1.258 },
285 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
286 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
287 1.442,1.400,1.354,1.319,1.299,1.272 },
288 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
289 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
290 1.456,1.412,1.364,1.328,1.307,1.282 }};
291
292 //data/corrections for T > Tlim
293
294 static const G4double hecorr[15] = {
295 120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
296 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
297 -22.30};
298
299 G4double sigma;
300 SetParticle(part);
301
302 G4double Z23 = G4Pow::GetInstance()->Z23(G4lrint(atomicNumber));
303
304 // correction if particle .ne. e-/e+
305 // compute equivalent kinetic energy
306 // lambda depends on p*beta ....
307
308 G4double eKineticEnergy = kinEnergy;
309
311 {
312 G4double TAU = kinEnergy/mass ;
313 G4double c = mass*TAU*(TAU+2.)/(CLHEP::electron_mass_c2*(TAU+1.)) ;
314 G4double w = c-2.;
315 G4double tau = 0.5*(w+std::sqrt(w*w+4.*c)) ;
316 eKineticEnergy = CLHEP::electron_mass_c2*tau ;
317 }
318
319 G4double eTotalEnergy = eKineticEnergy + CLHEP::electron_mass_c2 ;
320 G4double beta2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
321 /(eTotalEnergy*eTotalEnergy);
322 G4double bg2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
324
325 static const G4double epsfactor = 2.*CLHEP::electron_mass_c2*
328 G4double eps = epsfactor*bg2/Z23;
329
330 if (eps<epsmin) sigma = 2.*eps*eps;
331 else if(eps<epsmax) sigma = G4Log(1.+2.*eps)-2.*eps/(1.+2.*eps);
332 else sigma = G4Log(2.*eps)-1.+1./eps;
333
334 sigma *= chargeSquare*atomicNumber*atomicNumber/(beta2*bg2);
335
336 // interpolate in AtomicNumber and beta2
337 G4double c1,c2,cc1;
338
339 // get bin number in Z
340 G4int iZ = 14;
341 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
342 while ((iZ>=0)&&(Zdat[iZ]>=atomicNumber)) { --iZ; }
343
344 iZ = std::min(std::max(iZ, 0), 13);
345
346 G4double ZZ1 = Zdat[iZ];
347 G4double ZZ2 = Zdat[iZ+1];
348 G4double ratZ = (atomicNumber-ZZ1)*(atomicNumber+ZZ1)/
349 ((ZZ2-ZZ1)*(ZZ2+ZZ1));
350
351 static const G4double Tlim = 10.*CLHEP::MeV;
352 static const G4double sigmafactor =
354 static const G4double beta2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
356 static const G4double bg2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
358
359 static const G4double sig0[15] = {
360 0.2672*CLHEP::barn, 0.5922*CLHEP::barn, 2.653*CLHEP::barn, 6.235*CLHEP::barn,
361 11.69*CLHEP::barn , 13.24*CLHEP::barn , 16.12*CLHEP::barn, 23.00*CLHEP::barn,
362 35.13*CLHEP::barn , 39.95*CLHEP::barn , 50.85*CLHEP::barn, 67.19*CLHEP::barn,
363 91.15*CLHEP::barn , 104.4*CLHEP::barn , 113.1*CLHEP::barn};
364
365 static const G4double Tdat[22] = {
366 100*CLHEP::eV, 200*CLHEP::eV, 400*CLHEP::eV, 700*CLHEP::eV,
369 100*CLHEP::keV, 200*CLHEP::keV, 400*CLHEP::keV, 700*CLHEP::keV,
371 10*CLHEP::MeV, 20*CLHEP::MeV};
372
373 if(eKineticEnergy <= Tlim)
374 {
375 // get bin number in T (beta2)
376 G4int iT = 21;
377 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
378 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
379
380 iT = std::min(std::max(iT, 0), 20);
381
382 // calculate betasquare values
383 G4double T = Tdat[iT];
385 G4double b2small = T*(E+CLHEP::electron_mass_c2)/(E*E);
386
387 T = Tdat[iT+1];
389 G4double b2big = T*(E+CLHEP::electron_mass_c2)/(E*E);
390 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
391
392 if (charge < 0.)
393 {
394 c1 = celectron[iZ][iT];
395 c2 = celectron[iZ+1][iT];
396 cc1 = c1+ratZ*(c2-c1);
397
398 c1 = celectron[iZ][iT+1];
399 c2 = celectron[iZ+1][iT+1];
400 }
401 else
402 {
403 c1 = cpositron[iZ][iT];
404 c2 = cpositron[iZ+1][iT];
405 cc1 = c1+ratZ*(c2-c1);
406
407 c1 = cpositron[iZ][iT+1];
408 c2 = cpositron[iZ+1][iT+1];
409 }
410 G4double cc2 = c1+ratZ*(c2-c1);
411 sigma *= sigmafactor/(cc1+ratb2*(cc2-cc1));
412 }
413 else
414 {
415 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
416 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
417 if((atomicNumber >= ZZ1) && (atomicNumber <= ZZ2))
418 sigma = c1+ratZ*(c2-c1) ;
419 else if(atomicNumber < ZZ1)
420 sigma = atomicNumber*atomicNumber*c1/(ZZ1*ZZ1);
421 else if(atomicNumber > ZZ2)
422 sigma = atomicNumber*atomicNumber*c2/(ZZ2*ZZ2);
423 }
424 // low energy correction based on theory
425 sigma *= (1.+0.30/(1.+std::sqrt(1000.*eKineticEnergy)));
426
427 return sigma;
428}
429
430//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
431
433{
435 firstStep = true;
436 insideskin = false;
437 fr = facrange;
439 smallstep = 1.e10;
441 tlimitmin = 10.*tlimitminfix;
442 rndmEngineMod = G4Random::getTheEngine();
443}
444
445//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
446
448 const G4Track& track,
449 G4double& currentMinimalStep)
450{
451 tPathLength = currentMinimalStep;
452 const G4DynamicParticle* dp = track.GetDynamicParticle();
453
455 G4StepStatus stepStatus = sp->GetStepStatus();
458 idx = couple->GetIndex();
465 /*
466 G4cout << "G4Urban::StepLimit tPathLength= " << tPathLength
467 << " range= " <<currentRange<< " lambda= "<<lambda0
468 <<G4endl;
469 */
470
471 // stop here if small step
473 latDisplasment = false;
474 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
475 }
476
477 // upper limit for the straight line distance the particle can travel
478 // for electrons and positrons
479 G4double distance = (mass < masslimite)
480 ? currentRange*msc[idx]->doverra
481 // for muons, hadrons
482 : currentRange*msc[idx]->doverrb;
483
484 presafety = sp->GetSafety();
485 /*
486 G4cout << "G4Urban::StepLimit tPathLength= "
487 <<tPathLength<<" safety= " << presafety
488 << " range= " <<currentRange<< " lambda= "<<lambda0
489 << " Alg: " << steppingAlgorithm <<G4endl;
490 */
491 // far from geometry boundary
492 if(distance < presafety)
493 {
494 latDisplasment = false;
495 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
496 }
497
499 // standard version
500 //
502 {
503 //compute geomlimit and presafety
505 /*
506 G4cout << "G4Urban::Distance to boundary geomlimit= "
507 <<geomlimit<<" safety= " << presafety<<G4endl;
508 */
509
510 // is it far from boundary ?
511 if(distance < presafety)
512 {
513 latDisplasment = false;
514 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
515 }
516
517 smallstep += 1.;
518 insideskin = false;
519
520 // initialisation at firs step and at the boundary
521 if(firstStep || (stepStatus == fGeomBoundary))
522 {
524 if(!firstStep) { smallstep = 1.; }
525
526 //stepmin ~ lambda_elastic
530 /*
531 G4cout << "rangeinit= " << rangeinit << " stepmin= " << stepmin
532 << " tlimitmin= " << tlimitmin << " geomlimit= "
533 << geomlimit <<G4endl;
534 */
535 // constraint from the geometry
536
537 if((geomlimit < geombig) && (geomlimit > geommin))
538 {
539 // geomlimit is a geometrical step length
540 // transform it to true path length (estimation)
541 if(lambda0 > geomlimit) {
543 }
544 tgeom = (stepStatus == fGeomBoundary)
546 }
547 else
548 {
549 tgeom = geombig;
550 }
551 }
552
553 //step limit
556
557 //lower limit for tlimit
559 /*
560 G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit
561 << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
562 */
563 // shortcut
564 if((tPathLength < tlimit) && (tPathLength < presafety) &&
566 {
567 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
568 }
569
570 // step reduction near to boundary
571 if(smallstep <= skin)
572 {
573 tlimit = stepmin;
574 insideskin = true;
575 }
576 else if(geomlimit < geombig)
577 {
578 if(geomlimit > skindepth)
579 {
581 }
582 else
583 {
584 insideskin = true;
586 }
587 }
588
590
591 // randomise if not 'small' step and step determined by msc
595 }
596 // for 'normal' simulation with or without magnetic field
597 // there no small step/single scattering at boundaries
598 else if(steppingAlgorithm == fUseSafety)
599 {
600 if(stepStatus != fGeomBoundary) {
601 presafety = ComputeSafety(sp->GetPosition(),tPathLength);
602 }
603 /*
604 G4cout << "presafety= " << presafety
605 << " firstStep= " << firstStep
606 << " stepStatus= " << stepStatus
607 << G4endl;
608 */
609 // is far from boundary
610 if(distance < presafety)
611 {
612 latDisplasment = false;
613 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
614 }
615
616 if(firstStep || (stepStatus == fGeomBoundary)) {
618 fr = facrange;
619 // 9.1 like stepping for e+/e- only (not for muons,hadrons)
620 if(mass < masslimite)
621 {
623 if(lambda0 > lambdalimit) {
624 fr *= (0.75+0.25*lambda0/lambdalimit);
625 }
626 }
627 //lower limit for tlimit
630 }
631
632 //step limit
635
636 //lower limit for tlimit
638
639 // randomise if step determined by msc
642 }
643 // new stepping mode UseSafetyPlus
645 {
646 if(stepStatus != fGeomBoundary) {
647 presafety = ComputeSafety(sp->GetPosition(),tPathLength);
648 }
649 /*
650 G4cout << "presafety= " << presafety
651 << " firstStep= " << firstStep
652 << " stepStatus= " << stepStatus
653 << G4endl;
654 */
655 // is far from boundary
656 if(distance < presafety)
657 {
658 latDisplasment = false;
659 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
660 }
661
662 if(firstStep || (stepStatus == fGeomBoundary)) {
664 fr = facrange;
666 if(mass < masslimite)
667 {
668 rangecut = msc[idx]->ecut;
669 if(lambda0 > lambdalimit) {
670 fr *= (0.84+0.16*lambda0/lambdalimit);
671 }
672 }
673 //lower limit for tlimit
676 }
677 //step limit
680
681 //lower limit for tlimit
683
684 // condition for tPathLength from drr and finalr
685 if(currentRange > finalr) {
687 finalr*(1.-drr)*(2.-finalr/currentRange);
689 }
690
691 // condition safety
692 if(currentRange > rangecut) {
693 if(firstStep) {
695 } else if(stepStatus != fGeomBoundary && presafety > stepmin) {
697 }
698 }
699
700 // randomise if step determined by msc
703 }
704
705 // simple step limitation
706 else
707 {
708 if (stepStatus == fGeomBoundary)
709 {
713 }
714 // randomise if step determined by msc
717 }
718 firstStep = false;
719 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
720}
721
722//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
723
725{
727 par1 = -1. ;
728 par2 = par3 = 0. ;
729
730 // this correction needed to run MSC with eIoni and eBrem inactivated
731 // and makes no harm for a normal run
733
734 // do the true -> geom transformation
736
737 // z = t for very small tPathLength
739
740 /*
741 G4cout << "ComputeGeomPathLength: tpl= " << tPathLength
742 << " R= " << currentRange << " L0= " << lambda0
743 << " E= " << currentKinEnergy << " "
744 << particle->GetParticleName() << G4endl;
745 */
747
748 if ((tau <= tausmall) || insideskin) {
750
751 } else if (tPathLength < currentRange*dtrl) {
752 if(tau < taulim) zPathLength = tPathLength*(1.-0.5*tau) ;
753 else zPathLength = lambda0*(1.-G4Exp(-tau));
754
755 } else if(currentKinEnergy < mass || tPathLength == currentRange) {
756 par1 = 1./currentRange ;
757 par2 = 1./(par1*lambda0) ;
758 par3 = 1.+par2 ;
760 zPathLength =
762 } else {
763 zPathLength = 1./(par1*par3);
764 }
765
766 } else {
770
771 par1 = (lambda0-lambda1)/(lambda0*tPathLength);
772 //G4cout << "par1= " << par1 << " L1= " << lambda1 << G4endl;
773 par2 = 1./(par1*lambda0);
774 par3 = 1.+par2 ;
775 zPathLength = (1.-G4Exp(par3*G4Log(lambda1/lambda0)))/(par1*par3);
776 }
777
779 //G4cout<< "zPathLength= "<< zPathLength<< " L0= " << lambda0 << G4endl;
780 return zPathLength;
781}
782
783//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
784
786{
787 // step defined other than transportation
788 if(geomStepLength == zPathLength) {
789 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
790 // << " step= " << geomStepLength << " *** " << G4endl;
791 return tPathLength;
792 }
793
794 zPathLength = geomStepLength;
795
796 // t = z for very small step
797 if(geomStepLength < tlimitminfix2) {
798 tPathLength = geomStepLength;
799
800 // recalculation
801 } else {
802
803 G4double tlength = geomStepLength;
804 if((geomStepLength > lambda0*tausmall) && !insideskin) {
805
806 if(par1 < 0.) {
807 tlength = -lambda0*G4Log(1.-geomStepLength/lambda0) ;
808 } else {
809 const G4double par4 = par1*par3;
810 if(par4*geomStepLength < 1.) {
811 tlength = (1.-G4Exp(G4Log(1.-par4*geomStepLength)/par3))/par1;
812 } else {
813 tlength = currentRange;
814 }
815 }
816
817 if(tlength < geomStepLength) { tlength = geomStepLength; }
818 else if(tlength > tPathLength) { tlength = tPathLength; }
819 }
820 tPathLength = tlength;
821 }
822 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
823 // << " step= " << geomStepLength << " &&& " << G4endl;
824
825 return tPathLength;
826}
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
829
832 G4double /*safety*/)
833{
834 fDisplacement.set(0.0,0.0,0.0);
835 G4double kinEnergy = currentKinEnergy;
838 } else if(tPathLength > currentRange*0.01) {
841 }
842
843 if((kinEnergy <= CLHEP::eV) || (tPathLength <= tlimitminfix) ||
845
846 G4double cth = SampleCosineTheta(tPathLength,kinEnergy);
847
848 // protection against 'bad' cth values
849 if(std::abs(cth) >= 1.0) { return fDisplacement; }
850
851 G4double sth = std::sqrt((1.0 - cth)*(1.0 + cth));
853 G4ThreeVector newDirection(sth*std::cos(phi),sth*std::sin(phi),cth);
854 newDirection.rotateUz(oldDirection);
855
857 /*
858 G4cout << "G4UrbanMscModel::SampleSecondaries: e(MeV)= " << kineticEnergy
859 << " sinTheta= " << sth << " safety(mm)= " << safety
860 << " trueStep(mm)= " << tPathLength
861 << " geomStep(mm)= " << zPathLength
862 << G4endl;
863 */
864
866 if(dispAlg96) { SampleDisplacement(sth, phi); }
867 else { SampleDisplacementNew(cth, phi); }
868 fDisplacement.rotateUz(oldDirection);
869 }
870 return fDisplacement;
871}
872
873//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
874
876 G4double kinEnergy)
877{
878 G4double cth = 1.0;
879 G4double tau = trueStepLength/lambda0;
880
881 G4double lambda1 = GetTransportMeanFreePath(particle, kinEnergy);
882 if(std::abs(lambda1 - lambda0) > lambda0*0.01 && lambda1 > 0.)
883 {
884 // mean tau value
885 tau = trueStepLength*G4Log(lambda0/lambda1)/(lambda0-lambda1);
886 }
887
888 currentTau = tau;
889 lambdaeff = trueStepLength/currentTau;
891
892 if (tau >= taubig) { cth = -1.+2.*rndmEngineMod->flat(); }
893 else if (tau >= tausmall) {
894 static const G4double numlim = 0.01;
895 static const G4double onethird = 1./3.;
896 G4double xmeanth, x2meanth;
897 if(tau < numlim) {
898 xmeanth = 1.0 - tau*(1.0 - 0.5*tau);
899 x2meanth= 1.0 - tau*(5.0 - 6.25*tau)*onethird;
900 } else {
901 xmeanth = G4Exp(-tau);
902 x2meanth = (1.+2.*G4Exp(-2.5*tau))*onethird;
903 }
904
905 // too large step of low-energy particle
906 G4double relloss = 1. - kinEnergy/currentKinEnergy;
907 static const G4double rellossmax= 0.50;
908 if(relloss > rellossmax) {
909 return SimpleScattering(xmeanth,x2meanth);
910 }
911 // is step extreme small ?
912 G4bool extremesmallstep = false;
914 G4double theta0;
915 if(trueStepLength > tsmall) {
916 theta0 = ComputeTheta0(trueStepLength,kinEnergy);
917 } else {
918 theta0 = std::sqrt(trueStepLength/tsmall)
919 *ComputeTheta0(tsmall,kinEnergy);
920 extremesmallstep = true;
921 }
922
923 static const G4double onesixth = 1./6.;
924 static const G4double theta0max = CLHEP::pi*onesixth;
925 //G4cout << "Theta0= " << theta0 << " theta0max= " << theta0max
926 // << " sqrt(tausmall)= " << sqrt(tausmall) << G4endl;
927
928 // protection for very small angles
929 G4double theta2 = theta0*theta0;
930
931 if(theta2 < tausmall) { return cth; }
932
933 if(theta0 > theta0max) {
934 return SimpleScattering(xmeanth,x2meanth);
935 }
936
937 G4double x = theta2*(1.0 - theta2/12.);
938 if(theta2 > numlim) {
939 G4double sth = 2*std::sin(0.5*theta0);
940 x = sth*sth;
941 }
942
943 // parameter for tail
944 G4double ltau= G4Log(tau);
945 G4double u = !extremesmallstep ? G4Exp(ltau*onesixth)
946 : G4Exp(G4Log(tsmall/lambda0)*onesixth);
947
949 G4double xsi = msc[idx]->coeffc1 +
950 u*(msc[idx]->coeffc2+msc[idx]->coeffc3*u)+msc[idx]->coeffc4*xx;
951
952 // tail should not be too big
953 xsi = std::max(xsi, 1.9);
954 /*
955 if(KineticEnergy > 20*MeV && xsi < 1.6) {
956 G4cout << "G4UrbanMscModel::SampleCosineTheta: E(GeV)= "
957 << KineticEnergy/GeV
958 << " !!** c= " << xsi
959 << " **!! length(mm)= " << trueStepLength << " Zeff= " << Zeff
960 << " " << couple->GetMaterial()->GetName()
961 << " tau= " << tau << G4endl;
962 }
963 */
964
965 G4double c = xsi;
966
967 if(std::abs(c-3.) < 0.001) { c = 3.001; }
968 else if(std::abs(c-2.) < 0.001) { c = 2.001; }
969
970 G4double c1 = c-1.;
971 G4double ea = G4Exp(-xsi);
972 G4double eaa = 1.-ea ;
973 G4double xmean1 = 1.-(1.-(1.+xsi)*ea)*x/eaa;
974 G4double x0 = 1. - xsi*x;
975
976 // G4cout << " xmean1= " << xmean1 << " xmeanth= " << xmeanth << G4endl;
977
978 if(xmean1 <= 0.999*xmeanth) {
979 return SimpleScattering(xmeanth,x2meanth);
980 }
981 //from continuity of derivatives
982 G4double b = 1.+(c-xsi)*x;
983
984 G4double b1 = b+1.;
985 G4double bx = c*x;
986
987 G4double eb1 = G4Exp(G4Log(b1)*c1);
988 G4double ebx = G4Exp(G4Log(bx)*c1);
989 G4double d = ebx/eb1;
990
991 G4double xmean2 = (x0 + d - (bx - b1*d)/(c-2.))/(1. - d);
992
993 G4double f1x0 = ea/eaa;
994 G4double f2x0 = c1/(c*(1. - d));
995 G4double prob = f2x0/(f1x0+f2x0);
996
997 G4double qprob = xmeanth/(prob*xmean1+(1.-prob)*xmean2);
998
999 // sampling of costheta
1000 //G4cout << "c= " << c << " qprob= " << qprob << " eb1= " << eb1
1001 // << " c1= " << c1 << " b1= " << b1 << " bx= " << bx << " eb1= " << eb1
1002 // << G4endl;
1004 if(rndmarray[0] < qprob)
1005 {
1006 G4double var = 0;
1007 if(rndmarray[1] < prob) {
1008 cth = 1.+G4Log(ea+rndmEngineMod->flat()*eaa)*x;
1009 } else {
1010 var = (1.0 - d)*rndmEngineMod->flat();
1011 if(var < numlim*d) {
1012 var /= (d*c1);
1013 cth = -1.0 + var*(1.0 - 0.5*var*c)*(2. + (c - xsi)*x);
1014 } else {
1015 cth = 1. + x*(c - xsi - c*G4Exp(-G4Log(var + d)/c1));
1016 }
1017 }
1018 } else {
1019 cth = -1.+2.*rndmarray[1];
1020 }
1021 }
1022 return cth;
1023}
1024
1025//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1026
1028 G4double kinEnergy)
1029{
1030 // for all particles take the width of the central part
1031 // from a parametrization similar to the Highland formula
1032 // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
1033 G4double invbetacp = (kinEnergy+mass)/(kinEnergy*(kinEnergy+2.*mass));
1034 if(currentKinEnergy != kinEnergy) {
1035 invbetacp = std::sqrt(invbetacp*(currentKinEnergy+mass)/
1037 }
1038 G4double y = trueStepLength/currentRadLength;
1039
1040 if(particle == positron)
1041 {
1042 G4double Zeff = msc[idx]->Zeff;
1043 static const G4double xl= 0.6;
1044 static const G4double xh= 0.9;
1045 static const G4double e = 113.0;
1046 G4double corr;
1047
1048 G4double tau = std::sqrt(currentKinEnergy*kinEnergy)/mass;
1049 G4double x = std::sqrt(tau*(tau+2.)/((tau+1.)*(tau+1.)));
1050 G4double a = 0.994-4.08e-3*Zeff;
1051 G4double b = 7.16+(52.6+365./Zeff)/Zeff;
1052 G4double c = 1.000-4.47e-3*Zeff;
1053 G4double d = 1.21e-3*Zeff;
1054 if(x < xl) {
1055 corr = a*(1.-G4Exp(-b*x));
1056 } else if(x > xh) {
1057 corr = c+d*G4Exp(e*(x-1.));
1058 } else {
1059 G4double yl = a*(1.-G4Exp(-b*xl));
1060 G4double yh = c+d*G4Exp(e*(xh-1.));
1061 G4double y0 = (yh-yl)/(xh-xl);
1062 G4double y1 = yl-y0*xl;
1063 corr = y0*x+y1;
1064 }
1065 //==================================================================
1066 y *= corr*(1.+Zeff*(1.84035e-4*Zeff-1.86427e-2)+0.41125);
1067 }
1068
1069 static const G4double c_highland = 13.6*CLHEP::MeV;
1070 G4double theta0 = c_highland*std::abs(charge)*std::sqrt(y)*invbetacp;
1071
1072 // correction factor from e- scattering data
1073 theta0 *= (msc[idx]->coeffth1+msc[idx]->coeffth2*G4Log(y));
1074 return theta0;
1075}
1076
1077//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1078
1080{
1081 // simple and fast sampling
1082 // based on single scattering results
1083 // u = r/rmax : mean value
1084
1086 if(rmax > 0.)
1087 {
1088 G4double r = 0.73*rmax;
1089
1090 // simple distribution for v=Phi-phi=psi ~exp(-beta*v)
1091 // beta determined from the requirement that distribution should give
1092 // the same mean value than that obtained from the ss simulation
1093
1094 static const G4double cbeta = 2.160;
1095 static const G4double cbeta1 = 1. - G4Exp(-cbeta*CLHEP::pi);
1097 G4double psi = -G4Log(1. - rndmarray[0]*cbeta1)/cbeta;
1098 G4double Phi = (rndmarray[1] < 0.5) ? phi+psi : phi-psi;
1099 fDisplacement.set(r*std::cos(Phi),r*std::sin(Phi),0.0);
1100 }
1101}
1102
1103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1104
1106{
1107 // best sampling based on single scattering results
1108 G4double rmax =
1110 G4double r(0.0);
1111 G4double u(0.0);
1112 static const G4double reps = 5.e-3;
1113
1114 if(rmax > 0.)
1115 {
1116 static const G4double umax = 0.855;
1117 static const G4double wlow = 0.750;
1118
1119 static const G4double ralpha = 6.83e+0;
1120 static const G4double ra1 =-4.16179e+1;
1121 static const G4double ra2 = 1.12548e+2;
1122 static const G4double ra3 =-8.66665e+1;
1123 static const G4double ralpha1 = 0.751*ralpha;
1124 static const G4double ralpha2 =ralpha-ralpha1;
1125 static const G4double rwa1 = G4Exp(ralpha1*reps);
1126 static const G4double rwa2 = G4Exp(ralpha1*umax)-rwa1;
1127 static const G4double rejamax = 1.16456;
1128
1129 static const G4double rbeta = 2.18e+1;
1130 static const G4double rb0 = 4.81382e+2;
1131 static const G4double rb1 =-1.12842e+4;
1132 static const G4double rb2 = 4.57745e+4;
1133 static const G4double rbeta1 = 0.732*rbeta;
1134 static const G4double rbeta2 = rbeta-rbeta1;
1135 static const G4double rwb1 = G4Exp(-rbeta1*umax);
1136 static const G4double rwb2 = rwb1-G4Exp(-rbeta1*(1.-reps));
1137 static const G4double rejbmax = 1.62651;
1138
1139 G4int count = 0;
1140 G4double uc,rej;
1141
1142 if(rndmEngineMod->flat() < wlow)
1143 {
1144 do {
1146 u = G4Log(rwa1+rwa2*rndmarray[0])/ralpha1;
1147 uc = umax-u;
1148 rej = G4Exp(-ralpha2*uc)*
1149 (1.+ralpha*uc+ra1*uc*uc+ra2*uc*uc*uc+ra3*uc*uc*uc*uc);
1150 } while (rejamax*rndmarray[1] > rej && ++count < 1000);
1151 }
1152 else
1153 {
1154 do {
1156 u = -G4Log(rwb1-rwb2*rndmarray[0])/rbeta1;
1157 uc = u-umax;
1158 rej = G4Exp(-rbeta2*uc)*
1159 (1.+rbeta*uc+rb0*uc*uc+rb1*uc*uc*uc+rb2*uc*uc*uc*uc);
1160 } while (rejbmax*rndmarray[1] > rej && ++count < 1000);
1161 }
1162 r = rmax*u;
1163 }
1164
1165 if(r > 0.)
1166 {
1167 // sample Phi using lateral correlation
1168 // and r/rmax - (Phi-phi) correlation
1169 // v = Phi-phi = acos(latcorr/(r*sth))
1170 // from SS simulation f(v)*g(v)
1171 // f(v) ~ exp(-a1*v) normalized distribution
1172 // g(v) rejection function (0 < g(v) <= 1)
1173 G4double v, rej;
1174
1175 static const G4double peps = 1.e-4;
1176 static const G4double Pi = CLHEP::pi;
1177 static const G4double palpha[10] = {2.300e+0,2.490e+0,2.610e+0,2.820e+0,2.710e+0,
1178 2.750e+0,2.910e+0,3.400e+0,4.150e+0,5.400e+0};
1179 static const G4double palpha1[10]= {4.600e-2,1.245e-1,2.610e-1,2.820e-1,2.710e-1,
1180 6.875e-1,1.019e+0,1.360e+0,1.660e+0,2.430e+0};
1181 static const G4double pejmax[10] = {3.513,1.968,1.479,1.239,1.116,
1182 1.081,1.064,1.073,1.103,1.158};
1183
1184 static const G4double pa1[10] = { 3.218e+0, 2.412e+0, 2.715e+0, 2.787e+0, 2.541e+0,
1185 2.508e+0, 2.600e+0, 3.231e+0, 4.588e+0, 6.584e+0};
1186 static const G4double pa2[10] = {-5.528e-1, 2.523e+0, 1.738e+0, 2.082e+0, 1.423e+0,
1187 4.682e-1,-6.883e-1,-2.147e+0,-5.127e+0,-1.054e+1};
1188 static const G4double pa3[10] = { 3.618e+0, 2.032e+0, 2.341e+0, 2.172e+0, 7.205e-1,
1189 4.655e-1, 6.318e-1, 1.255e+0, 2.425e+0, 4.938e+0};
1190 static const G4double pa4[10] = { 2.437e+0, 9.450e-1, 4.349e-1, 2.221e-1, 1.130e-1,
1191 5.405e-2, 2.245e-2, 7.370e-3, 1.456e-3, 1.508e-4};
1192 static const G4double pw1[10] = {G4Exp(-palpha1[0]*peps),G4Exp(-palpha1[1]*peps),
1193 G4Exp(-palpha1[2]*peps),G4Exp(-palpha1[3]*peps),
1194 G4Exp(-palpha1[4]*peps),G4Exp(-palpha1[5]*peps),
1195 G4Exp(-palpha1[6]*peps),G4Exp(-palpha1[7]*peps),
1196 G4Exp(-palpha1[8]*peps),G4Exp(-palpha1[9]*peps)};
1197 static const G4double pw2[10] = {pw1[0]-G4Exp(-palpha1[0]*(Pi-peps)),
1198 pw1[1]-G4Exp(-palpha1[1]*(Pi-peps)),
1199 pw1[2]-G4Exp(-palpha1[2]*(Pi-peps)),
1200 pw1[3]-G4Exp(-palpha1[3]*(Pi-peps)),
1201 pw1[4]-G4Exp(-palpha1[4]*(Pi-peps)),
1202 pw1[5]-G4Exp(-palpha1[5]*(Pi-peps)),
1203 pw1[6]-G4Exp(-palpha1[6]*(Pi-peps)),
1204 pw1[7]-G4Exp(-palpha1[7]*(Pi-peps)),
1205 pw1[8]-G4Exp(-palpha1[8]*(Pi-peps)),
1206 pw1[9]-G4Exp(-palpha1[9]*(Pi-peps))};
1207
1208 G4int iphi = (G4int)(u*10.);
1209 if(iphi < 0) { iphi = 0; }
1210 else if(iphi > 9) { iphi = 9; }
1211 G4int count = 0;
1212
1213 do {
1215 v = -G4Log(pw1[iphi]-pw2[iphi]*rndmarray[0])/palpha1[iphi];
1216 rej = (G4Exp(-palpha[iphi]*v)*
1217 (1+pa1[iphi]*v+pa2[iphi]*v*v+pa3[iphi]*v*v*v)+pa4[iphi])/
1218 G4Exp(-pw1[iphi]*v);
1219 }
1220 // Loop checking, 5-March-2018, Vladimir Ivanchenko
1221 while (pejmax[iphi]*rndmarray[1] > rej && ++count < 1000);
1222
1223 G4double Phi = (rndmEngineMod->flat() < 0.5) ? phi+v : phi-v;
1224 fDisplacement.set(r*std::cos(Phi),r*std::sin(Phi),0.0);
1225 }
1226}
1227
1228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1229
1231{
1232 // it is assumed, that for the second run only addition
1233 // of a new G4MaterialCutsCouple is possible
1234 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
1235 size_t numOfCouples = theCoupleTable->GetTableSize();
1236 if(numOfCouples != msc.size()) { msc.resize(numOfCouples, nullptr); }
1237
1238 for(size_t j=0; j<numOfCouples; ++j) {
1239 auto aCouple = theCoupleTable->GetMaterialCutsCouple(j);
1240
1241 // cut may be changed before runs
1242 G4double cut = aCouple->GetProductionCuts()->GetProductionCut(1);
1243 if(nullptr != msc[j]) {
1244 msc[j]->ecut = cut;
1245 continue;
1246 }
1247 // new couple
1248 msc[j] = new mscData();
1249 msc[j]->ecut = cut;
1250 G4double Zeff = aCouple->GetMaterial()->GetIonisation()->GetZeffective();
1251 msc[j]->Zeff = Zeff;
1252 msc[j]->sqrtZ = std::sqrt(Zeff);
1253 G4double lnZ = G4Log(Zeff);
1254 // correction in theta0 formula
1255 G4double w = G4Exp(lnZ/6.);
1256 G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
1257 msc[j]->coeffth1 = facz*(1. - 8.7780e-2/Zeff);
1258 msc[j]->coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
1259
1260 // tail parameters
1261 G4double Z13 = w*w;
1262 msc[j]->coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
1263 msc[j]->coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
1264 msc[j]->coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
1265 msc[j]->coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
1266
1267 msc[j]->Z23 = Z13*Z13;
1268
1269 msc[j]->stepmina = 27.725/(1.+0.203*Zeff);
1270 msc[j]->stepminb = 6.152/(1.+0.111*Zeff);
1271
1272 // 21.07.2020
1273 msc[j]->doverra = 9.6280e-1 - 8.4848e-2*msc[j]->sqrtZ + 4.3769e-3*Zeff;
1274
1275 // 06.10.2020
1276 // msc[j]->doverra = 7.7024e-1 - 6.7878e-2*msc[j]->sqrtZ + 3.5015e-3*Zeff;
1277
1278 msc[j]->doverrb = 1.15 - 9.76e-4*Zeff;
1279 }
1280}
1281
1282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const G4double eps
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
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 &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
G4double GetLogKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4bool LateralDisplacementAlg96() const
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
const G4Material * GetMaterial() const
G4double GetRadlen() const
Definition: G4Material.hh:216
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static G4Positron * Positron()
Definition: G4Positron.cc:93
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
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
void SampleDisplacementNew(G4double sinTheta, G4double phi)
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double ComputeTrueStepLength(G4double geomStepLength) override
G4double currentKinEnergy
void SampleDisplacement(G4double sinTheta, G4double phi)
G4double rndmarray[2]
G4double Randomizetlimit()
const G4ParticleDefinition * positron
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
void StartTracking(G4Track *) override
~G4UrbanMscModel() override
G4double currentRadLength
const G4MaterialCutsCouple * couple
G4UrbanMscModel(const G4String &nam="UrbanMsc")
G4double currentLogKinEnergy
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeTlimitmin()
static std::vector< mscData * > msc
G4double ComputeGeomPathLength(G4double truePathLength) override
void SetParticle(const G4ParticleDefinition *)
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
G4ParticleChangeForMSC * fParticleChange
const G4ParticleDefinition * particle
CLHEP::HepRandomEngine * rndmEngineMod
G4LossTableManager * theManager
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
G4double ComputeStepmin()
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
G4double dtrl
Definition: G4VMscModel.hh:200
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:159
G4double facrange
Definition: G4VMscModel.hh:196
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:290
G4double skin
Definition: G4VMscModel.hh:199
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:319
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:78
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:224
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:189
G4double lambdalimit
Definition: G4VMscModel.hh:201
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:209
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:280
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:272
G4double facsafety
Definition: G4VMscModel.hh:198
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:208
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:116
G4double facgeom
Definition: G4VMscModel.hh:197
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double electron_mass_c2
static constexpr double proton_mass_c2
static constexpr double um
Definition: SystemOfUnits.h:94
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double keV
static constexpr double twopi
Definition: SystemOfUnits.h:56
static constexpr double MeV
static constexpr double Bohr_radius
static constexpr double classic_electr_radius
static constexpr double eV
static constexpr double hbarc
static constexpr double pi
Definition: SystemOfUnits.h:55
static constexpr double nm
Definition: SystemOfUnits.h:93
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define LOG_EKIN_MIN
Definition: templates.hh:98
int G4lrint(double ad)
Definition: templates.hh:134