Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UrbanMscModel.hh
Go to the documentation of this file.
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: $
27 // GEANT4 tag $Name: $
28 //
29 // -------------------------------------------------------------------
30 //
31 //
32 // GEANT4 Class header file
33 //
34 //
35 // File name: G4UrbanMscModel
36 //
37 // Author: Laszlo Urban
38 //
39 // Creation date: 19.02.2013
40 //
41 // Created from G4UrbanMscModel96
42 //
43 // New parametrization for theta0
44 // Correction for very small step length
45 //
46 // Class Description:
47 //
48 // Implementation of the model of multiple scattering based on
49 // H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4UrbanMscModel_h
55 #define G4UrbanMscModel_h 1
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 
61 #include "G4VMscModel.hh"
62 #include "G4MscStepLimitType.hh"
63 #include "G4Log.hh"
64 #include "G4Exp.hh"
65 
67 class G4SafetyHelper;
68 class G4LossTableManager;
69 
70 static const G4double c_highland = 13.6*CLHEP::MeV ;
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 {
76 
77 public:
78 
79  G4UrbanMscModel(const G4String& nam = "UrbanMsc");
80 
81  virtual ~G4UrbanMscModel();
82 
83  void Initialise(const G4ParticleDefinition*, const G4DataVector&);
84 
85  void StartTracking(G4Track*);
86 
88  G4double KineticEnergy,
89  G4double AtomicNumber,
90  G4double AtomicWeight=0.,
91  G4double cut =0.,
92  G4double emax=DBL_MAX);
93 
95 
97  G4double& currentMinimalStep);
98 
99  G4double ComputeGeomPathLength(G4double truePathLength);
100 
101  G4double ComputeTrueStepLength(G4double geomStepLength);
102 
103  inline G4double ComputeTheta0(G4double truePathLength,
104  G4double KineticEnergy);
105 
106 private:
107 
108  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
109 
110  inline void SetParticle(const G4ParticleDefinition*);
111 
112  inline void UpdateCache();
113 
114  inline G4double SampleDisplacement();
115 
116  inline G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
117 
118  inline G4double LatCorrelation();
119 
120  // hide assignment operator
121  G4UrbanMscModel & operator=(const G4UrbanMscModel &right);
123 
124  const G4ParticleDefinition* particle;
125  G4ParticleChangeForMSC* fParticleChange;
126 
127  const G4MaterialCutsCouple* couple;
128  G4LossTableManager* theManager;
129 
130  G4double mass;
131  G4double charge,ChargeSquare;
132  G4double masslimite,lambdalimit,fr;
133 
134  G4double taubig;
135  G4double tausmall;
136  G4double taulim;
137  G4double currentTau;
138  G4double tlimit;
139  G4double tlimitmin;
140  G4double tlimitminfix,tlimitminfix2;
141  G4double tgeom;
142 
143  G4double geombig;
144  G4double geommin;
145  G4double geomlimit;
146  G4double skindepth;
147  G4double smallstep;
148 
149  G4double presafety;
150 
151  G4double lambda0;
152  G4double lambdaeff;
153  G4double tPathLength;
154  G4double zPathLength;
155  G4double par1,par2,par3;
156 
157  G4double stepmin;
158 
159  G4double currentKinEnergy;
160  G4double currentRange;
161  G4double rangeinit;
162  G4double currentRadLength;
163 
164  G4double theta0max,rellossmax;
165  G4double third;
166 
167  G4int currentMaterialIndex;
168 
169  G4double y;
170  G4double Zold;
171  G4double Zeff,Z2,Z23,lnZ;
172  G4double coeffth1,coeffth2;
173  G4double coeffc1,coeffc2,coeffc3,coeffc4;
174 
175  G4bool firstStep;
176  G4bool inside;
177  G4bool insideskin;
178 };
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182 
183 inline
184 void G4UrbanMscModel::SetParticle(const G4ParticleDefinition* p)
185 {
186  if (p != particle) {
187  particle = p;
188  mass = p->GetPDGMass();
189  charge = p->GetPDGCharge()/CLHEP::eplus;
190  ChargeSquare = charge*charge;
191  }
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 
196 inline
197 void G4UrbanMscModel::UpdateCache()
198 {
199  lnZ = G4Log(Zeff);
200  // correction in theta0 formula
201  G4double w = G4Exp(lnZ/6.);
202  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
203  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
204  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
205 
206  // tail parameters
207  G4double Z13 = std::exp(lnZ/3.);
208  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
209  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
210  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
211  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
212 
213  Z2 = Zeff*Zeff;
214  Z23 = Z13*Z13;
215 
216  Zold = Zeff;
217 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 
221 inline
223  G4double KineticEnergy)
224 {
225  // for all particles take the width of the central part
226  // from a parametrization similar to the Highland formula
227  // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
228  G4double invbetacp = sqrt((currentKinEnergy+mass)*(KineticEnergy+mass)/
229  (currentKinEnergy*(currentKinEnergy+2.*mass)*
230  KineticEnergy*(KineticEnergy+2.*mass)));
231  y = trueStepLength/currentRadLength;
232  G4double theta0 = c_highland*std::abs(charge)*sqrt(y)*invbetacp;
233  y = G4Log(y);
234  // correction factor from e- scattering data
235  theta0 *= (coeffth1+coeffth2*y);
236  return theta0;
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240 
241 inline
242 G4double G4UrbanMscModel::SimpleScattering(G4double xmeanth, G4double x2meanth)
243 {
244  // 'large angle scattering'
245  // 2 model functions with correct xmean and x2mean
246  G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
247  G4double prob = (a+2.)*xmeanth/a;
248 
249  // sampling
250  G4double cth = 1.;
251  if(G4UniformRand() < prob) {
252  cth = -1.+2.*G4Exp(G4Log(G4UniformRand())/(a+1.));
253  } else {
254  cth = -1.+2.*G4UniformRand();
255  }
256  return cth;
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260 
261 inline
262 G4double G4UrbanMscModel::SampleDisplacement()
263 {
264  G4double r = 0.0;
265  if ((currentTau >= tausmall) && !insideskin) {
266  G4double rmax = sqrt((tPathLength-zPathLength)*(tPathLength+zPathLength));
267  r = rmax*G4Exp(G4Log(G4UniformRand())*third);
268  }
269  return r;
270 }
271 
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
273 
274 inline
275 G4double G4UrbanMscModel::LatCorrelation()
276 {
277  static const G4double kappa = 2.5;
278  static const G4double kappami1 = 1.5;
279 
280  G4double latcorr = 0.;
281  if((currentTau >= tausmall) && !insideskin)
282  {
283  if(currentTau < taulim)
284  latcorr = lambdaeff*kappa*currentTau*currentTau*
285  (1.-(kappa+1.)*currentTau*third)*third;
286  else
287  {
288  G4double etau = 0.;
289  if(currentTau < taubig) { etau = G4Exp(-currentTau); }
290  latcorr = -kappa*currentTau;
291  latcorr = G4Exp(latcorr)/kappami1;
292  latcorr += 1.-kappa*etau/kappami1 ;
293  latcorr *= 2.*lambdaeff*third;
294  }
295  }
296  return latcorr;
297 }
298 
299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
300 
301 #endif
302 
G4UrbanMscModel(const G4String &nam="UrbanMsc")
virtual ~G4UrbanMscModel()
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
bool G4bool
Definition: G4Types.hh:79
G4double ComputeGeomPathLength(G4double truePathLength)
G4double ComputeTrueStepLength(G4double geomStepLength)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
void StartTracking(G4Track *)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX)
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)