source: trunk/source/processes/electromagnetic/standard/src/G4UrbanMscModel90.cc @ 1340

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

update ti head

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