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