Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UrbanMscModel92.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: G4UrbanMscModel92.hh 66592 2012-12-23 09:34:55Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4UrbanMscModel92
35 //
36 // Author: Laszlo Urban
37 //
38 // Creation date: 06.03.2008
39 //
40 // Modifications:
41 //
42 // 28-10-2009 V.Ivanchenko moved G4UrbanMscModel to G4UrbanMscModel92,
43 // now it is a frozen version of the Urban model corresponding
44 // to g4 9.2
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 G4UrbanMscModel92_h
55 #define G4UrbanMscModel92_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  G4UrbanMscModel92(const G4String& nam = "UrbanMsc92");
76 
77  virtual ~G4UrbanMscModel92();
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  G4UrbanMscModel92 & operator=(const G4UrbanMscModel92 &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;
170  G4double scr1ini,scr2ini,scr1,scr2;
171 
172  G4bool firstStep;
173  G4bool inside;
174  G4bool insideskin;
175 };
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 
180 inline
181 void G4UrbanMscModel92::SetParticle(const G4ParticleDefinition* p)
182 {
183  if (p != particle) {
184  particle = p;
185  mass = p->GetPDGMass();
186  charge = p->GetPDGCharge()/CLHEP::eplus;
187  ChargeSquare = charge*charge;
188  }
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
193 inline
194 void G4UrbanMscModel92::UpdateCache()
195 {
196  lnZ = std::log(Zeff);
197  coeffth1 = 0.885+lnZ*(0.104-0.0170*lnZ);
198  coeffth2 = 0.028+lnZ*(0.012-0.00125*lnZ);
199  coeffc1 = 2.134-lnZ*(0.1045-0.00602*lnZ);
200  coeffc2 = 0.001126-lnZ*(0.0001089+0.0000247*lnZ);
201  Z2 = Zeff*Zeff;
202  Z23 = std::exp(2.*lnZ/3.);
203  scr1 = scr1ini*Z23;
204  scr2 = scr2ini*Z2*ChargeSquare;
205  // lastMaterial = couple->GetMaterial();
206  Zold = Zeff;
207 }
208 
209 #endif
210