source: trunk/source/processes/electromagnetic/standard/src/G4UrbanMscModel.cc @ 1228

Last change on this file since 1228 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

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