source: trunk/source/processes/electromagnetic/standard/src/G4MscModel71.cc @ 991

Last change on this file since 991 was 991, checked in by garnier, 15 years ago

update

File size: 25.9 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: G4MscModel71.cc,v 1.6 2008/03/13 17:20:07 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:   G4MscModel71
35//
36// Author:      Laszlo Urban
37//
38// Creation date: 03.03.2001
39//
40// Modifications:
41//
42// 27-03-03 Move model part from G4MultipleScattering (V.Ivanchenko)
43// 23-05-03 important change in angle distribution for muons/hadrons
44//          the central part now is similar to the Highland parametrization +
45//          minor correction in angle sampling algorithm (for all particles)
46//          (L.Urban)
47// 30-05-03 misprint in SampleCosineTheta corrected(L.Urban)
48// 27-03-03 Rename (V.Ivanchenko)
49// 05-08-03 angle distribution has been modified (L.Urban)
50// 06-11-03 precision problems solved for high energy (PeV) particles
51//          change in the tail of the angular distribution
52//          highKinEnergy is set to 100 PeV (L.Urban)
53//
54// 10-11-03 highKinEnergy is set back to 100 TeV, some tail tuning +
55//          cleaning (L.Urban)
56// 26-11-03 correction in TrueStepLength :
57//          trueLength <= currentRange (L.Urban)
58// 01-03-04 signature changed in SampleCosineTheta,
59//          energy dependence calculations has been simplified,
60// 11-03-04 corrections in GeomPathLength,TrueStepLength,
61//          SampleCosineTheta
62// 23-04-04 true -> geom and geom -> true transformation has been
63//          rewritten, changes in the angular distribution (L.Urban)
64// 19-07-04 correction in SampleCosineTheta in order to avoid
65//          num. precision problems at high energy/small step(L.Urban)
66// 17-08-04 changes in the angle distribution (slightly modified
67//          Highland formula for the width of the central part,
68//          changes in the numerical values of some other parameters)
69//          ---> approximately step independent distribution (L.Urban)
70// 21-09-04 change in the tail of the angular distribution (L.Urban)
71//
72// 03-11-04 precision problem for very high energy ions and small stepsize
73//          solved in SampleCosineTheta (L.Urban).
74// 15-04-05 optimize internal interface - add SampleSecondaries method (V.Ivanchenko)
75// 03-10-05 Model is freezed with the name McsModel71 (V.Ivanchenko)
76// 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
77//
78
79// Class Description:
80//
81// Implementation of the model of multiple scattering based on
82// H.W.Lewis Phys Rev 78 (1950) 526 and others
83
84// -------------------------------------------------------------------
85//
86
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
91#include "G4MscModel71.hh"
92#include "Randomize.hh"
93#include "G4Electron.hh"
94#include "G4LossTableManager.hh"
95#include "G4PhysicsTable.hh"
96#include "G4ParticleChangeForMSC.hh"
97#include "G4TransportationManager.hh"
98#include "G4Navigator.hh"
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
102using namespace std;
103
104G4MscModel71::G4MscModel71(G4double& m_dtrl, G4double& m_NuclCorrPar,
105                           G4double& m_FactPar, G4double& m_factail,
106                       G4bool& m_samplez, const G4String& nam)
107  : G4VEmModel(nam),
108  taubig(8.0),
109  tausmall(1.e-20),
110  taulim(1.e-6),
111  dtrl(m_dtrl),
112  NuclCorrPar (m_NuclCorrPar),
113  FactPar(m_FactPar),
114  factail(m_factail),
115  samplez(m_samplez),
116  isInitialized(false)
117{
118  stepmin       = 1.e-6*mm;
119  currentRange  = 0. ;
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
124G4MscModel71::~G4MscModel71()
125{}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
129void G4MscModel71::Initialise(const G4ParticleDefinition* p,
130                              const G4DataVector&)
131{
132  if(isInitialized) return;
133  // set values of some data members
134  sigmafactor = twopi*classic_electr_radius*classic_electr_radius;
135  particle = p;
136  mass = particle->GetPDGMass();
137  charge = particle->GetPDGCharge()/eplus;
138  b = 1. ;
139  xsi = 3.00 ;
140
141  if(pParticleChange)
142    fParticleChange = reinterpret_cast<G4ParticleChangeForMSC*>(pParticleChange);
143  else
144    fParticleChange = new G4ParticleChangeForMSC();
145
146  navigator = G4TransportationManager::GetTransportationManager()
147    ->GetNavigatorForTracking();
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151
152G4double G4MscModel71::ComputeCrossSectionPerAtom(
153                             const G4ParticleDefinition* part,
154                                   G4double KineticEnergy,
155                                   G4double AtomicNumber,
156                                   G4double AtomicWeight, 
157                                   G4double,
158                                   G4double)
159{
160  const G4double epsfactor = 2.*electron_mass_c2*electron_mass_c2*
161                             Bohr_radius*Bohr_radius/(hbarc*hbarc);
162  const G4double epsmin = 1.e-4 , epsmax = 1.e10;
163
164  const G4double Zdat[15] = { 4., 6.,13.,20.,26.,29.,32.,38.,47.,
165                             50.,56.,64.,74.,79.,82. };
166
167  const G4double Tdat[23] = {0.0001*MeV,0.0002*MeV,0.0004*MeV,0.0007*MeV,
168                             0.001*MeV,0.002*MeV,0.004*MeV,0.007*MeV,
169                             0.01*MeV,0.02*MeV,0.04*MeV,0.07*MeV,
170                             0.1*MeV,0.2*MeV,0.4*MeV,0.7*MeV,
171                             1.*MeV,2.*MeV,4.*MeV,7.*MeV,10.*MeV,20.*MeV,
172                             10000.0*MeV};
173
174 // corr. factors for e-/e+ lambda
175          G4double celectron[15][23] =
176          {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
177            1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
178            1.112,1.108,1.100,1.093,1.089,1.087,0.7235     },
179           {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
180            1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
181            1.109,1.105,1.097,1.090,1.086,1.082,0.7925     },
182           {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
183            1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
184            1.131,1.124,1.113,1.104,1.099,1.098,0.9147     },
185           {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
186            1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
187            1.112,1.105,1.096,1.089,1.085,1.098,0.9700     },
188           {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
189            1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
190            1.073,1.070,1.064,1.059,1.056,1.056,1.0022     },
191           {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
192            1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
193            1.074,1.070,1.063,1.059,1.056,1.052,1.0158     },
194           {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
195            1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
196            1.068,1.064,1.059,1.054,1.051,1.050,1.0284     },
197           {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
198            1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
199            1.039,1.037,1.034,1.031,1.030,1.036,1.0515     },
200           {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
201            1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
202            1.031,1.028,1.024,1.022,1.021,1.024,1.0834     },
203           {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
204            1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
205            1.020,1.017,1.015,1.013,1.013,1.020,1.0937     },
206           {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
207            1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
208            0.995,0.993,0.993,0.993,0.993,1.011,1.1140     },
209           {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
210            1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
211            0.974,0.972,0.973,0.974,0.975,0.987,1.1410     },
212           {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
213            1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
214            0.950,0.947,0.949,0.952,0.954,0.963,1.1750     },
215           {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
216            1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
217            0.941,0.938,0.940,0.944,0.946,0.954,1.1922     },
218           {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
219            1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
220            0.933,0.930,0.933,0.936,0.939,0.949,1.2026     }};
221           G4double cpositron[15][23] = {
222           {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
223            1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
224            1.131,1.126,1.117,1.108,1.103,1.100,0.7235     },
225           {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
226            1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
227            1.138,1.132,1.122,1.113,1.108,1.102,0.7925     },
228           {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
229            1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
230            1.203,1.190,1.173,1.159,1.151,1.145,0.9147     },
231           {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
232            1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
233            1.225,1.210,1.191,1.175,1.166,1.174,0.9700     },
234           {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
235            1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
236            1.217,1.203,1.184,1.169,1.160,1.151,1.0022     },
237           {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
238            1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
239            1.237,1.222,1.201,1.184,1.174,1.159,1.0158     },
240           {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
241            1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
242            1.252,1.234,1.212,1.194,1.183,1.170,1.0284     },
243           {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
244            2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
245            1.254,1.237,1.214,1.195,1.185,1.179,1.0515     },
246           {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
247            2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
248            1.312,1.288,1.258,1.235,1.221,1.205,1.0834     },
249           {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
250            2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
251            1.320,1.294,1.264,1.240,1.226,1.214,1.0937     },
252           {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
253            2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
254            1.328,1.302,1.270,1.245,1.231,1.233,1.1140     },
255           {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
256            2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
257            1.361,1.330,1.294,1.267,1.251,1.239,1.1410     },
258           {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
259            3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
260            1.409,1.372,1.330,1.298,1.280,1.258,1.1750     },
261           {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
262            3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
263            1.442,1.400,1.354,1.319,1.299,1.272,1.1922     },
264           {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
265            3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
266            1.456,1.412,1.364,1.328,1.307,1.282,1.2026     }};
267
268  G4double sigma;
269  if (part != particle ) {
270    particle = part;
271    mass = particle->GetPDGMass();
272    charge = particle->GetPDGCharge()/eplus;
273  }
274
275  G4double Z23 = 2.*log(AtomicNumber)/3.; Z23 = exp(Z23);
276
277  // correction if particle .ne. e-/e+
278  // compute equivalent kinetic energy
279  // lambda depends on p*beta ....
280
281  G4double eKineticEnergy = KineticEnergy;
282
283  if((particle->GetParticleName() != "e-") &&
284     (particle->GetParticleName() != "e+") )
285  {
286     G4double TAU = KineticEnergy/mass ;
287     G4double c = mass*TAU*(TAU+2.)/(electron_mass_c2*(TAU+1.)) ;
288     G4double w = c-2. ;
289     G4double tau = 0.5*(w+sqrt(w*w+4.*c)) ;
290     eKineticEnergy = electron_mass_c2*tau ;
291  }
292
293  G4double ChargeSquare = charge*charge;
294
295  G4double eTotalEnergy = eKineticEnergy + electron_mass_c2 ;
296  G4double beta2 = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
297                                 /(eTotalEnergy*eTotalEnergy);
298  G4double bg2   = eKineticEnergy*(eTotalEnergy+electron_mass_c2)
299                                 /(electron_mass_c2*electron_mass_c2);
300
301  G4double eps = epsfactor*bg2/Z23;
302
303  if     (eps<epsmin)  sigma = 2.*eps*eps;
304  else if(eps<epsmax)  sigma = log(1.+2.*eps)-2.*eps/(1.+2.*eps);
305  else                 sigma = log(2.*eps)-1.+1./eps;
306
307  sigma *= ChargeSquare*AtomicNumber*AtomicNumber/(beta2*bg2);
308
309  // nuclear size effect correction for high energy
310  // ( a simple approximation at present)
311  G4double corrnuclsize,a,w1,w2,w;
312
313  G4double x0 = 1. - NuclCorrPar*mass/(KineticEnergy*
314               exp(log(AtomicWeight/(g/mole))/3.));
315  if ( x0 < -1. || eKineticEnergy  <= 10.*MeV)
316      {
317        x0 = -1.;
318        corrnuclsize = 1.;
319      }
320  else
321      {
322        a = 1.+1./eps;
323        if (eps > epsmax) w1=log(2.*eps)+1./eps-3./(8.*eps*eps);
324        else              w1=log((a+1.)/(a-1.))-2./(a+1.);
325        w = 1./((1.-x0)*eps);
326        if (w < epsmin)   w2=-log(w)-1.+2.*w-1.5*w*w;
327        else              w2 = log((a-x0)/(a-1.))-(1.-x0)/(a-x0);
328        corrnuclsize = w1/w2;
329        corrnuclsize = exp(-FactPar*mass/KineticEnergy)*
330                      (corrnuclsize-1.)+1.;
331      }
332
333  // interpolate in AtomicNumber and beta2
334  // get bin number in Z
335  G4int iZ = 14;
336  while ((iZ>=0)&&(Zdat[iZ]>=AtomicNumber)) iZ -= 1;
337  if (iZ==14)                               iZ = 13;
338  if (iZ==-1)                               iZ = 0 ;
339
340  G4double Z1 = Zdat[iZ];
341  G4double Z2 = Zdat[iZ+1];
342  G4double ratZ = (AtomicNumber-Z1)/(Z2-Z1);
343
344  // get bin number in T (beta2)
345  G4int iT = 22;
346  while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
347  if(iT==22)                                  iT = 21;
348  if(iT==-1)                                  iT = 0 ;
349
350  //  calculate betasquare values
351  G4double T = Tdat[iT],   E = T + electron_mass_c2;
352  G4double b2small = T*(E+electron_mass_c2)/(E*E);
353  T = Tdat[iT+1]; E = T + electron_mass_c2;
354  G4double b2big = T*(E+electron_mass_c2)/(E*E);
355  G4double ratb2 = (beta2-b2small)/(b2big-b2small);
356  G4double c1,c2,cc1,cc2,corr;
357  if (charge < 0.)
358    {
359       c1 = celectron[iZ][iT];
360       c2 = celectron[iZ+1][iT];
361       cc1 = c1+ratZ*(c2-c1);
362
363       c1 = celectron[iZ][iT+1];
364       c2 = celectron[iZ+1][iT+1];
365       cc2 = c1+ratZ*(c2-c1);
366
367       corr = cc1+ratb2*(cc2-cc1);
368       sigma /= corr;
369    }
370
371  if (charge > 0.)
372    {
373       c1 = cpositron[iZ][iT];
374       c2 = cpositron[iZ+1][iT];
375       cc1 = c1+ratZ*(c2-c1);
376
377       c1 = cpositron[iZ][iT+1];
378       c2 = cpositron[iZ+1][iT+1];
379       cc2 = c1+ratZ*(c2-c1);
380
381       corr = cc1+ratb2*(cc2-cc1);
382       sigma /= corr;
383    }
384
385  sigma *= sigmafactor;
386
387  //  nucl. size correction for particles other than e+/e- only at present !!!!
388  if((particle->GetParticleName() != "e-") &&
389     (particle->GetParticleName() != "e+")   )
390     sigma /= corrnuclsize;
391  //  G4cout << "e= " << KineticEnergy << " sigma= " << sigma << G4endl;
392
393  return sigma;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
398G4double G4MscModel71::GeomPathLength(
399                          G4PhysicsTable* theLambdaTable,
400                    const G4MaterialCutsCouple* couple,
401                    const G4ParticleDefinition* theParticle,
402                          G4double& T0,
403                          G4double lambda,
404                          G4double range,
405                          G4double truePathLength)
406{
407  //  do the true -> geom transformation
408  const G4double  ztmax = 101./103. ;
409  if (theParticle != particle ) {
410    particle = theParticle;
411    mass = particle->GetPDGMass();
412    charge = particle->GetPDGCharge()/eplus;
413  }
414  currentKinEnergy = T0;
415  currentRange = range ;
416  currentRadLength = couple->GetMaterial()->GetRadlen();
417
418  lambda0 = lambda;
419  par1 = -1. ; 
420  par2 = par3 = 0. ; 
421  tPathLength = truePathLength;
422
423  // this correction needed to run MSC with eIoni and eBrem inactivated
424  // and makes no harm for a normal run
425  if(tPathLength > range)
426    tPathLength = range ;
427
428  G4double tau   = tPathLength/lambda0 ;
429
430  if (tau <= tausmall) return tPathLength;
431
432  G4double zmean = tPathLength;
433  if (tPathLength < range*dtrl) {
434    zmean = lambda0*(1.-exp(-tau));
435    if(tau < taulim) zmean = tPathLength*(1.-0.5*tPathLength/lambda0) ;
436  } else if(T0 < mass) {
437    par1 = 1./range ;
438    par2 = 1./(par1*lambda0) ;
439    par3 = 1.+par2 ;
440    zmean = (1.-exp(par3*log(1.-tPathLength/range)))/(par1*par3) ;
441  } else {
442    G4LossTableManager* theManager = G4LossTableManager::Instance();
443    G4double T1 = theManager->GetEnergy(particle,range-tPathLength,couple);
444    G4double lambda1 ;
445    if (theLambdaTable) {
446      G4bool bb;
447      lambda1 = ((*theLambdaTable)[couple->GetIndex()])->GetValue(T1,bb);
448    } else {
449      lambda1 = CrossSection(couple,particle,T1,0.0,1.0);
450    }
451    lambda1 = 1.0/lambda1;
452    par1 = (lambda0-lambda1)/(lambda0*tPathLength) ;
453    par2 = 1./(par1*lambda0) ;
454    par3 = 1.+par2 ;
455    zmean = (1.-exp(par3*log(lambda1/lambda0)))/(par1*par3) ;
456  }
457
458  //  sample z
459  G4double zPathLength = zmean ;
460  G4double zt = zmean/tPathLength ;
461  if (tPathLength >= stepmin && samplez && zt > 0.5 && zt < ztmax)
462  {
463    G4double cz = 0.5*(3.*zt-1.)/(1.-zt) ;
464    G4double cz1 = 1.+cz ;
465    G4double u0 = cz/cz1 ;
466    G4double u,grej ;
467    do {
468        u = exp(log(G4UniformRand())/cz1) ;
469        grej = exp(cz*log(u/u0))*(1.-u)/(1.-u0) ;
470      } while (grej < G4UniformRand()) ;
471   zPathLength = tPathLength*u ;
472  }
473
474  return zPathLength ;
475}
476
477//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
478
479G4double G4MscModel71::TrueStepLength(G4double geomStepLength)
480{
481  G4double trueLength = geomStepLength;
482  trueLength = geomStepLength;
483  if(geomStepLength > lambda0*tausmall)
484  {
485    if(par1 <  0.)
486      trueLength = -lambda0*log(1.-geomStepLength/lambda0) ;
487    else 
488    {
489      if(par1*par3*geomStepLength < 1.)
490        trueLength = (1.-exp(log(1.-par1*par3*geomStepLength)/par3))/par1 ;
491      else 
492        trueLength = currentRange ;
493    } 
494  }
495
496  if(trueLength > tPathLength) trueLength = tPathLength;
497  if(trueLength > currentRange) trueLength = currentRange ;
498  if(trueLength < geomStepLength) trueLength = geomStepLength;
499
500  return trueLength;
501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
505void G4MscModel71::SampleSecondaries(std::vector<G4DynamicParticle*>*,
506                                     const G4MaterialCutsCouple*,
507                                     const G4DynamicParticle* dynParticle,
508                                     G4double truestep,
509                                     G4double safety)
510{
511  G4double kineticEnergy = dynParticle->GetKineticEnergy();
512  if(kineticEnergy <= 0.0) return;
513
514  G4double cth  = SampleCosineTheta(truestep,kineticEnergy);
515  G4double sth  = sqrt((1.0 - cth)*(1.0 + cth));
516  G4double phi  = twopi*G4UniformRand();
517  G4double dirx = sth*cos(phi);
518  G4double diry = sth*sin(phi);
519
520  G4ThreeVector oldDirection = dynParticle->GetMomentumDirection();
521  G4ThreeVector newDirection(dirx,diry,cth);
522  newDirection.rotateUz(oldDirection);
523  fParticleChange->ProposeMomentumDirection(newDirection);
524
525  /*
526    const G4ParticleDefinition* pd = dynParticle->GetDefinition();
527    G4cout << "G4MscModel71: Sample secondary; E(MeV)= " << kineticEnergy/MeV
528           << " MeV; step(mm)= " << truestep/mm
529           << ", safety(mm)= " <<  safety/mm << " " << pd->GetParticleName()
530           << G4endl;
531  */
532
533  if (latDisplasment && safety > 0.0) {
534
535    G4double r = SampleDisplacement();
536    if (r > safety) r = safety;
537
538    // sample direction of lateral displacement
539    G4double phi  = twopi*G4UniformRand();
540    G4double dirx = std::cos(phi);
541    G4double diry = std::sin(phi);
542
543    G4ThreeVector newPosition(dirx,diry,0.0);
544    newPosition.rotateUz(oldDirection);
545
546    // compute new endpoint of the Step
547    newPosition *= r;
548    newPosition += *(fParticleChange->GetProposedPosition());
549
550    navigator->LocateGlobalPointWithinVolume(newPosition);
551
552    fParticleChange->ProposePosition(newPosition);
553  }
554}
555
556//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
557
558G4double G4MscModel71::SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
559{
560  G4double cth = 1. ;
561  G4double tau = trueStepLength/lambda0 ;
562
563  if(trueStepLength >= currentRange*dtrl) {
564    if(par1*trueStepLength < 1.)
565      tau = -par2*log(1.-par1*trueStepLength) ;
566    else
567      tau = taubig ;
568  }
569  currentTau = tau ;
570
571  if(trueStepLength < stepmin)
572    cth = exp(-tau) ;
573  else
574  {
575    if (tau >= taubig) cth = -1.+2.*G4UniformRand();
576    else if (tau >= tausmall)
577    {
578        G4double a ;
579
580        // for all particles take the width of the central part
581        //  from a  parametrization similar to the Highland formula
582        // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
583        // here : theta0 = 13.6*MeV*Q*(t/X0)**0.555/(beta*cp)
584        const G4double c_highland = 13.6*MeV, corr_highland=0.555 ;
585        G4double Q = std::abs(charge) ;
586        G4double xx0 = trueStepLength/currentRadLength;
587        G4double betacp = sqrt(currentKinEnergy*(currentKinEnergy+2.*mass)*
588                               KineticEnergy*(KineticEnergy+2.*mass)/
589                               ((currentKinEnergy+mass)*(KineticEnergy+mass))) ;
590        G4double theta0 = c_highland*Q*exp(corr_highland*log(xx0))/betacp ;
591
592        if(theta0 > taulim) a = 0.5/(1.-cos(theta0)) ;
593        else                a = 1.0/(theta0*theta0) ;
594
595        G4double xmeanth = exp(-tau);
596        G4double xmeanth1 = 1.-xmeanth ;
597        if(currentTau < taulim) xmeanth1 = tau ;         
598
599        const G4double x1fac1 = exp(-xsi) ;
600        const G4double x1fac2 = (1.-(1.+xsi)*x1fac1)/(1.-x1fac1) ;
601        const G4double x1fac3 = 1.3      ; 
602
603        G4double ea,eaa,xmean1 ;
604        G4double c = 2.,b1 = 2., bx = 2.,
605                 eb1 = b1, ebx = b1, xmean2 = 0. ;
606        G4double prob = 1., qprob ;             
607        G4double x0 = 1.-xsi/a;
608        G4double oneminusx0=xsi/a ;
609        G4double oneplusx0=2.+xsi/a ;
610
611        G4double f1x0=1., f2x0=1. ;
612        const G4double tau0 = 0.10 ;
613        if(tau > tau0)
614        {
615         // 1 model function
616          a = 1./xmeanth1 ;
617          ea = exp(-2.*a) ;
618          eaa= 1.-ea ;
619          xmean1 = 1.-1./a+2.*ea/eaa ;
620          prob = 1. ;
621          qprob = 1. ;
622        }
623        else if (x0 <= -1.)
624        {
625          // 2 model fuctions only
626          // in order to have xmean1 > xmeanth -> qprob < 1
627          x0 = -1.;
628
629          if( a < 1./xmeanth1)
630            a = 1./xmeanth1 ;
631         
632          oneminusx0 = 1.-x0 ;
633          oneplusx0 =  1.+x0 ;
634          ea = exp(-a*oneminusx0);
635          eaa = 1.-ea ;
636          xmean1 = 1.-1./a+oneminusx0*ea/eaa ;
637          qprob = xmeanth/xmean1 ;
638
639        }
640        else
641        {
642          // 3 model fuctions
643          // in order to have xmean1 > xmeanth
644          if((1.-x1fac2/a) < xmeanth)
645          {
646            a = x1fac3*x1fac2/xmeanth1 ;
647            x0 = 1.-xsi/a ;
648            oneminusx0=xsi/a ;
649            oneplusx0=2.-xsi/a ;
650          }
651
652          ea = x1fac1 ;
653          eaa = 1.-ea ;
654          xmean1 = 1.-x1fac2/a ;
655
656          const G4double fctail = factail*1.0 ;
657          c = 2.+fctail*tau ;
658          G4double c1 = c-1. ;
659          G4double c2 = c-2. ;
660          if(c2 == 0.) c2 = fctail*tausmall ;           
661
662          b = 1.+(c-xsi)/a ;
663
664          b1 = b+1. ;
665          bx = c/a ;
666          eb1=exp((c1)*log(b1)) ;
667          ebx=exp((c1)*log(bx)) ;
668          xmean2 = (x0*eb1+ebx-(eb1*bx-b1*ebx)/c2)/(eb1-ebx) ;
669
670          f1x0 = a*ea/eaa ;
671          f2x0 = c1*eb1*ebx/(eb1-ebx)/
672                          exp(c*log(bx)) ;
673          // from continuity at x=x0
674          prob = f2x0/(f1x0+f2x0) ;
675          // from xmean = xmeanth
676          qprob = (f1x0+f2x0)*xmeanth/(f2x0*xmean1+f1x0*xmean2) ;
677
678        }
679
680        // sampling of costheta
681        if (G4UniformRand() < qprob)
682        {
683          if (G4UniformRand() < prob)
684             cth = 1.+log(ea+G4UniformRand()*eaa)/a ;
685          else
686             cth = b-b1*bx/exp(log(ebx-G4UniformRand()*(ebx-eb1))/(c-1.)) ;
687        }
688        else
689        {
690          cth = -1.+2.*G4UniformRand();
691        }
692    }
693  } 
694
695  return cth ;
696}
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
699G4double G4MscModel71::SampleDisplacement()
700{
701  const G4double kappa = 2.5;
702  const G4double kappapl1 = kappa+1.;
703  const G4double kappami1 = kappa-1.;
704  G4double rmean = 0.0;
705  if (currentTau >= tausmall) {
706    if (currentTau < taulim) {
707      rmean = kappa*currentTau*currentTau*currentTau*(1.-kappapl1*currentTau*0.25)/6. ;
708
709    } else {
710      G4double etau = 0.0;
711      if (currentTau<taubig) etau = exp(-currentTau);
712      rmean = -kappa*currentTau;
713      rmean = -exp(rmean)/(kappa*kappami1);
714      rmean += currentTau-kappapl1/kappa+kappa*etau/kappami1;
715    }
716    if (rmean>0.) rmean = 2.*lambda0*sqrt(rmean/3.0);
717    else          rmean = 0.;
718 }
719  return rmean;
720}
721
722//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
723
724
725
726
Note: See TracBrowser for help on using the repository browser.