Geant4_10
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 {
74 
75 public:
76 
77  G4UrbanMscModel(const G4String& nam = "UrbanMsc");
78 
79  virtual ~G4UrbanMscModel();
80 
81  void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82 
83  void StartTracking(G4Track*);
84 
86  G4double KineticEnergy,
87  G4double AtomicNumber,
88  G4double AtomicWeight=0.,
89  G4double cut =0.,
90  G4double emax=DBL_MAX);
91 
93 
95  G4double& currentMinimalStep);
96 
97  G4double ComputeGeomPathLength(G4double truePathLength);
98 
99  G4double ComputeTrueStepLength(G4double geomStepLength);
100 
101  G4double ComputeTheta0(G4double truePathLength,
102  G4double KineticEnergy);
103 
104 private:
105 
106  G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
107 
108  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
109 
110  G4double SampleDisplacement();
111 
112  G4double LatCorrelation();
113 
114  inline void SetParticle(const G4ParticleDefinition*);
115 
116  inline void UpdateCache();
117 
118  // hide assignment operator
119  G4UrbanMscModel & operator=(const G4UrbanMscModel &right);
121 
122  const G4ParticleDefinition* particle;
123  G4ParticleChangeForMSC* fParticleChange;
124 
125  const G4MaterialCutsCouple* couple;
126  G4LossTableManager* theManager;
127 
128  G4double mass;
129  G4double charge,ChargeSquare;
130  G4double masslimite,lambdalimit,fr;
131 
132  G4double taubig;
133  G4double tausmall;
134  G4double taulim;
135  G4double currentTau;
136  G4double tlimit;
137  G4double tlimitmin;
138  G4double tlimitminfix,tlimitminfix2;
139  G4double tgeom;
140 
141  G4double geombig;
142  G4double geommin;
143  G4double geomlimit;
144  G4double skindepth;
145  G4double smallstep;
146 
147  G4double presafety;
148 
149  G4double lambda0;
150  G4double lambdaeff;
151  G4double tPathLength;
152  G4double zPathLength;
153  G4double par1,par2,par3;
154 
155  G4double stepmin;
156 
157  G4double currentKinEnergy;
158  G4double currentRange;
159  G4double rangeinit;
160  G4double currentRadLength;
161 
162  G4double theta0max,rellossmax;
163  G4double third;
164 
165  G4int currentMaterialIndex;
166 
167  G4double y;
168  G4double Zold;
169  G4double Zeff,Z2,Z23,lnZ;
170  G4double coeffth1,coeffth2;
171  G4double coeffc1,coeffc2,coeffc3,coeffc4;
172 
173  G4bool firstStep;
174  G4bool inside;
175  G4bool insideskin;
176 };
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
181 inline
182 void G4UrbanMscModel::SetParticle(const G4ParticleDefinition* p)
183 {
184  if (p != particle) {
185  particle = p;
186  mass = p->GetPDGMass();
187  charge = p->GetPDGCharge()/CLHEP::eplus;
188  ChargeSquare = charge*charge;
189  }
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
194 inline
195 void G4UrbanMscModel::UpdateCache()
196 {
197  lnZ = G4Log(Zeff);
198  // correction in theta0 formula
199  G4double w = G4Exp(lnZ/6.);
200  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
201  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
202  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
203 
204  // tail parameters
205  G4double Z13 = std::exp(lnZ/3.);
206  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
207  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
208  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
209  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
210 
211  Z2 = Zeff*Zeff;
212  Z23 = Z13*Z13;
213 
214  Zold = Zeff;
215 }
216 
217 #endif
218 
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
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)