source: trunk/source/processes/electromagnetic/standard/src/G4UrbanMscModel2.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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