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

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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