Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UrbanMscModel96.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: G4UrbanMscModel96
36 //
37 // Author: Laszlo Urban
38 //
39 // Creation date: 26.09.2012
40 //
41 // Created from G4UrbanMscModel95
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 G4UrbanMscModel96_h
55 #define G4UrbanMscModel96_h 1
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 
61 #include "G4VMscModel.hh"
62 #include "G4MscStepLimitType.hh"
63 
65 class G4SafetyHelper;
66 class G4LossTableManager;
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69 
71 {
72 
73 public:
74 
75  G4UrbanMscModel96(const G4String& nam = "UrbanMsc96");
76 
77  virtual ~G4UrbanMscModel96();
78 
79  void Initialise(const G4ParticleDefinition*, const G4DataVector&);
80 
81  void StartTracking(G4Track*);
82 
84  G4double KineticEnergy,
85  G4double AtomicNumber,
86  G4double AtomicWeight=0.,
87  G4double cut =0.,
88  G4double emax=DBL_MAX);
89 
91 
93  G4double& currentMinimalStep);
94 
95  G4double ComputeGeomPathLength(G4double truePathLength);
96 
97  G4double ComputeTrueStepLength(G4double geomStepLength);
98 
99  G4double ComputeTheta0(G4double truePathLength,
100  G4double KineticEnergy);
101 
102 private:
103 
104  G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
105 
106  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
107 
108  G4double SampleDisplacement();
109 
110  G4double LatCorrelation();
111 
112  inline void SetParticle(const G4ParticleDefinition*);
113 
114  inline void UpdateCache();
115 
116  // hide assignment operator
117  G4UrbanMscModel96 & operator=(const G4UrbanMscModel96 &right);
119 
120  const G4ParticleDefinition* particle;
121  G4ParticleChangeForMSC* fParticleChange;
122 
123  const G4MaterialCutsCouple* couple;
124  G4LossTableManager* theManager;
125 
126  G4double mass;
127  G4double charge,ChargeSquare;
128  G4double masslimite,lambdalimit,fr;
129 
130  G4double taubig;
131  G4double tausmall;
132  G4double taulim;
133  G4double currentTau;
134  G4double tlimit;
135  G4double tlimitmin;
136  G4double tlimitminfix;
137  G4double tgeom;
138 
139  G4double geombig;
140  G4double geommin;
141  G4double geomlimit;
142  G4double skindepth;
143  G4double smallstep;
144 
145  G4double presafety;
146 
147  G4double lambda0;
148  G4double lambdaeff;
149  G4double tPathLength;
150  G4double zPathLength;
151  G4double par1,par2,par3;
152 
153  G4double stepmin;
154 
155  G4double currentKinEnergy;
156  G4double currentRange;
157  G4double rangeinit;
158  G4double currentRadLength;
159 
160  G4double theta0max,rellossmax;
161  G4double third;
162 
163  G4int currentMaterialIndex;
164 
165  G4double y;
166  G4double Zold;
167  G4double Zeff,Z2,Z23,lnZ;
168  G4double coeffth1,coeffth2;
169  G4double coeffc1,coeffc2,coeffc3,coeffc4;
170 
171  G4bool firstStep;
172  G4bool inside;
173  G4bool insideskin;
174 };
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178 
179 inline
180 void G4UrbanMscModel96::SetParticle(const G4ParticleDefinition* p)
181 {
182  if (p != particle) {
183  particle = p;
184  mass = p->GetPDGMass();
185  charge = p->GetPDGCharge()/CLHEP::eplus;
186  ChargeSquare = charge*charge;
187  }
188 }
189 
190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191 
192 inline
193 void G4UrbanMscModel96::UpdateCache()
194 {
195  lnZ = std::log(Zeff);
196  // correction in theta0 formula
197  G4double facz = 0.74845+0.13354*exp(log(Zeff)/6.);
198  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
199  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
200 
201  // tail parameters
202  G4double Z13 = std::exp(lnZ/3.);
203  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
204  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
205  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
206  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
207 
208  Z2 = Zeff*Zeff;
209  Z23 = Z13*Z13;
210 
211  Zold = Zeff;
212 }
213 
214 #endif
215