Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UrbanMscModel93.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//
28// -------------------------------------------------------------------
29//
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4UrbanMscModel93
35//
36// Author: Laszlo Urban
37//
38// Creation date: 06.03.2008
39//
40// Modifications:
41// 23-04-2009 L.Urban updated parameterization in UpdateCache method
42// 28-10-2009 V.Ivanchenko moved G4UrbanMscModel2 to G4UrbanMscModel93,
43// now it is a frozen version of the Urban model corresponding
44// to g4 9.3
45//
46//
47// Class Description:
48//
49// Implementation of the model of multiple scattering based on
50// H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
51
52// -------------------------------------------------------------------
53//
54
55#ifndef G4UrbanMscModel93_h
56#define G4UrbanMscModel93_h 1
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
61
62#include "G4VMscModel.hh"
63#include "G4MscStepLimitType.hh"
64
66class G4SafetyHelper;
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
72{
73
74public:
75
76 G4UrbanMscModel93(const G4String& nam = "UrbanMsc93");
77
78 virtual ~G4UrbanMscModel93();
79
80 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
81
83
85 G4double KineticEnergy,
86 G4double AtomicNumber,
87 G4double AtomicWeight=0.,
88 G4double cut =0.,
89 G4double emax=DBL_MAX);
90
92
94 G4double& currentMinimalStep);
95
97
99
100 G4double ComputeTheta0(G4double truePathLength,
101 G4double KineticEnergy);
102
103private:
104
105 G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
106
107 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
108
109 G4double SampleDisplacement();
110
111 G4double LatCorrelation();
112
113 inline void SetParticle(const G4ParticleDefinition*);
114
115 inline void UpdateCache();
116
117 // hide assignment operator
118 G4UrbanMscModel93 & operator=(const G4UrbanMscModel93 &right);
120
121 const G4ParticleDefinition* particle;
122 G4ParticleChangeForMSC* fParticleChange;
123
124 const G4MaterialCutsCouple* couple;
125 G4LossTableManager* theManager;
126
127 G4double mass;
128 G4double charge,ChargeSquare;
129 G4double masslimite,lambdalimit,fr;
130
131 G4double taubig;
132 G4double tausmall;
133 G4double taulim;
134 G4double currentTau;
135 G4double tlimit;
136 G4double tlimitmin;
137 G4double tlimitminfix;
138 G4double tgeom;
139
140 G4double geombig;
141 G4double geommin;
142 G4double geomlimit;
143 G4double skindepth;
144 G4double smallstep;
145
146 G4double presafety;
147
148 G4double lambda0;
149 G4double lambdaeff;
150 G4double tPathLength;
151 G4double zPathLength;
152 G4double par1,par2,par3;
153
154 G4double stepmin;
155
156 G4double currentKinEnergy;
157 G4double currentRange;
158 G4double rangeinit;
159 G4double currentRadLength;
160
161 G4double theta0max,rellossmax;
162 G4double third;
163
164 G4int currentMaterialIndex;
165
166 G4double y;
167 G4double Zold;
168 G4double Zeff,Z2,Z23,lnZ;
169 G4double coeffth1,coeffth2;
170 G4double coeffc1,coeffc2;
171 G4double scr1ini,scr2ini,scr1,scr2;
172
173 G4bool firstStep;
174 G4bool inside;
175 G4bool insideskin;
176};
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
181inline
182void G4UrbanMscModel93::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
194inline
195void G4UrbanMscModel93::UpdateCache()
196{
197 lnZ = std::log(Zeff);
198 //new correction in theta0 formula
199 coeffth1 = (1.-8.7780e-2/Zeff)*(0.87+0.03*lnZ);
200 coeffth2 = (4.0780e-2+1.7315e-4*Zeff)*(0.87+0.03*lnZ);
201 // tail parameters
202 G4double lnZ1 = std::log(Zeff+1.);
203 coeffc1 = 2.943-0.197*lnZ1;
204 coeffc2 = 0.0987-0.0143*lnZ1;
205 // for single scattering
206 Z2 = Zeff*Zeff;
207 Z23 = std::exp(2.*lnZ/3.);
208 scr1 = scr1ini*Z23;
209 scr2 = scr2ini*Z2*ChargeSquare;
210
211 Zold = Zeff;
212}
213
214#endif
215
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetPDGCharge() const
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
void StartTracking(G4Track *)
G4double ComputeGeomPathLength(G4double truePathLength)
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)
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double ComputeTrueStepLength(G4double geomStepLength)
#define DBL_MAX
Definition: templates.hh:83