Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UrbanMscModel95.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: G4UrbanMscModel95
36//
37// Author: Laszlo Urban
38//
39// Creation date: 20.03.2011
40//
41// Created from G4UrbanMscModel93
42//
43// Class Description:
44//
45// Implementation of the model of multiple scattering based on
46// H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
47
48// -------------------------------------------------------------------
49//
50
51#ifndef G4UrbanMscModel95_h
52#define G4UrbanMscModel95_h 1
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
57
58#include "G4VMscModel.hh"
59#include "G4MscStepLimitType.hh"
60
62class G4SafetyHelper;
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
68{
69
70public:
71
72 G4UrbanMscModel95(const G4String& nam = "UrbanMsc95");
73
74 virtual ~G4UrbanMscModel95();
75
76 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
77
79
81 G4double KineticEnergy,
82 G4double AtomicNumber,
83 G4double AtomicWeight=0.,
84 G4double cut =0.,
85 G4double emax=DBL_MAX);
86
88
90 G4double& currentMinimalStep);
91
93
95
96 G4double ComputeTheta0(G4double truePathLength,
97 G4double KineticEnergy);
98
99private:
100
101 G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
102
103 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
104
105 G4double SampleDisplacement();
106
107 G4double LatCorrelation();
108
109 inline void SetParticle(const G4ParticleDefinition*);
110
111 inline void UpdateCache();
112
113 // hide assignment operator
114 G4UrbanMscModel95 & operator=(const G4UrbanMscModel95 &right);
116
117 const G4ParticleDefinition* particle;
118 G4ParticleChangeForMSC* fParticleChange;
119
120 const G4MaterialCutsCouple* couple;
121 G4LossTableManager* theManager;
122
123 G4double mass;
124 G4double charge,ChargeSquare;
125 G4double masslimite,lambdalimit,fr;
126
127 G4double taubig;
128 G4double tausmall;
129 G4double taulim;
130 G4double currentTau;
131 G4double tlimit;
132 G4double tlimitmin;
133 G4double tlimitminfix;
134 G4double tgeom;
135
136 G4double geombig;
137 G4double geommin;
138 G4double geomlimit;
139 G4double skindepth;
140 G4double smallstep;
141
142 G4double presafety;
143
144 G4double lambda0;
145 G4double lambdaeff;
146 G4double tPathLength;
147 G4double zPathLength;
148 G4double par1,par2,par3;
149
150 G4double stepmin;
151
152 G4double currentKinEnergy;
153 G4double currentRange;
154 G4double rangeinit;
155 G4double currentRadLength;
156
157 G4double theta0max,rellossmax;
158 G4double third;
159
160 G4int currentMaterialIndex;
161
162 G4double y;
163 G4double Zold;
164 G4double Zeff,Z2,Z23,lnZ;
165 G4double coeffth1,coeffth2;
166 G4double coeffc1,coeffc2,coeffc3,coeffc4;
167 G4double scr1ini,scr2ini,scr1,scr2;
168
169 G4bool firstStep;
170 G4bool inside;
171 G4bool insideskin;
172};
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176
177inline
178void G4UrbanMscModel95::SetParticle(const G4ParticleDefinition* p)
179{
180 if (p != particle) {
181 particle = p;
182 mass = p->GetPDGMass();
183 charge = p->GetPDGCharge()/CLHEP::eplus;
184 ChargeSquare = charge*charge;
185 }
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189
190inline
191void G4UrbanMscModel95::UpdateCache()
192{
193 lnZ = std::log(Zeff);
194 // correction in theta0 formula
195 coeffth1 = (1. - 8.7780e-2/Zeff)*(0.87 + 0.03*lnZ);
196 coeffth2 = (4.0780e-2 + 1.7315e-4*Zeff)*(0.87 + 0.03*lnZ);
197
198 // tail parameters
199 G4double Z13 = std::exp(lnZ/3.);
200 coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
201 coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
202 coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
203 coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
204
205 // for single scattering
206 Z2 = Zeff*Zeff;
207 Z23 = Z13*Z13;
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
G4double ComputeTrueStepLength(G4double geomStepLength)
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
void StartTracking(G4Track *)
G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX)
G4double ComputeGeomPathLength(G4double truePathLength)
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
#define DBL_MAX
Definition: templates.hh:83