Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27// -------------------------------------------------------------------
28//
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4UrbanMscModel
34//
35// Author: Laszlo Urban
36//
37// Creation date: 19.02.2013
38//
39// Created from G4UrbanMscModel96
40//
41// New parametrization for theta0
42// Correction for very small step length
43//
44// Class Description:
45//
46// Implementation of the model of multiple scattering based on
47// H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
48
49// -------------------------------------------------------------------
50//
51
52#ifndef G4UrbanMscModel_h
53#define G4UrbanMscModel_h 1
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
58
59#include "G4VMscModel.hh"
60#include "G4MscStepLimitType.hh"
61#include "G4Log.hh"
62#include "G4Exp.hh"
63
65class G4SafetyHelper;
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
70{
71
72public:
73
74 explicit G4UrbanMscModel(const G4String& nam = "UrbanMsc");
75
76 ~G4UrbanMscModel() override;
77
79 const G4DataVector&) override;
80
81 void StartTracking(G4Track*) override;
82
85 G4double KineticEnergy,
86 G4double AtomicNumber,
87 G4double AtomicWeight=0.,
88 G4double cut =0.,
89 G4double emax=DBL_MAX) override;
90
92 G4double safety) override;
93
95 G4double& currentMinimalStep) override;
96
97 G4double ComputeGeomPathLength(G4double truePathLength) override;
98
99 G4double ComputeTrueStepLength(G4double geomStepLength) override;
100
101 G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy);
102
103 // hide assignment operator
104 G4UrbanMscModel & operator=(const G4UrbanMscModel &right) = delete;
106
107private:
108
109 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
110
111 void SampleDisplacement(G4double sinTheta, G4double phi);
112
113 void SampleDisplacementNew(G4double sinTheta, G4double phi);
114
115 void InitialiseModelCache();
116
117 inline void SetParticle(const G4ParticleDefinition*);
118
119 inline G4double Randomizetlimit();
120
121 inline G4double SimpleScattering();
122
123 inline G4double ComputeStepmin();
124
125 inline G4double ComputeTlimitmin();
126
127 CLHEP::HepRandomEngine* rndmEngineMod;
128
129 const G4ParticleDefinition* particle;
130 const G4ParticleDefinition* positron;
131 G4ParticleChangeForMSC* fParticleChange;
132 const G4MaterialCutsCouple* couple;
133
134 G4double mass;
135 G4double charge,chargeSquare;
136 G4double masslimite,fr;
137
138 G4double taubig;
139 G4double tausmall;
140 G4double taulim;
141 G4double currentTau;
142 G4double tlimit;
143 G4double tlimitmin;
144 G4double tlimitminfix,tlimitminfix2;
145 G4double tgeom;
146
147 G4double geombig;
148 G4double geommin;
149 G4double geomlimit;
150 G4double skindepth;
151 G4double smallstep;
152
153 G4double presafety;
154
155 G4double lambda0;
156 G4double lambdaeff;
157 G4double tPathLength;
158 G4double zPathLength;
159 G4double par1,par2,par3;
160
161 G4double stepmin;
162
163 G4double currentKinEnergy;
164 G4double currentLogKinEnergy;
165 G4double currentRange;
166 G4double rangeinit;
167 G4double currentRadLength;
168
169 G4double drr,finalr;
170
171 G4double tlow;
172 G4double invmev;
173 G4double xmeanth = 0.0;
174 G4double x2meanth = 1./3.;
175 G4double rndmarray[2];
176
177 struct mscData {
178 G4double Z23, sqrtZ;
179 G4double coeffth1, coeffth2;
180 G4double coeffc1, coeffc2, coeffc3, coeffc4;
181 G4double stepmina, stepminb;
182 G4double doverra, doverrb;
183 G4double posa, posb, posc, posd, pose;
184 };
185 static std::vector<mscData*> msc;
186
187 // index of G4MaterialCutsCouple
188 G4int idx;
189
190 G4bool firstStep;
191 G4bool insideskin;
192
193 G4bool latDisplasmentbackup;
194 G4bool dispAlg96;
195 G4bool fPosiCorrection = true;
196 G4bool isFirstInstance = false;
197};
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
201
202inline
203void G4UrbanMscModel::SetParticle(const G4ParticleDefinition* p)
204{
205 if (p != particle) {
206 particle = p;
207 mass = p->GetPDGMass();
208 charge = p->GetPDGCharge()/CLHEP::eplus;
209 chargeSquare = charge*charge;
210 }
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214
215inline G4double G4UrbanMscModel::Randomizetlimit()
216{
217 G4double res = tlimitmin;
218 if(tlimit > tlimitmin)
219 {
220 res = G4RandGauss::shoot(rndmEngineMod,tlimit,0.1*(tlimit-tlimitmin));
221 res = std::max(res, tlimitmin);
222 }
223 return res;
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227
228inline G4double G4UrbanMscModel::SimpleScattering()
229{
230 // 'large angle scattering'
231 // 2 model functions with correct xmean and x2mean
232 G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
233 G4double prob = (a+2.)*xmeanth/a;
234
235 // sampling
236 rndmEngineMod->flatArray(2, rndmarray);
237 return (rndmarray[1] < prob) ?
238 -1.+2.*G4Exp(G4Log(rndmarray[0])/(a+1.)) : -1.+2.*rndmarray[0];
239}
240
241inline G4double G4UrbanMscModel::ComputeStepmin()
242{
243 // define stepmin using estimation of the ratio
244 // of lambda_elastic/lambda_transport
245 G4double rat = currentKinEnergy*invmev;
246 return lambda0*1.e-3/(2.e-3+rat*(msc[idx]->stepmina+msc[idx]->stepminb*rat));
247}
248
249inline G4double G4UrbanMscModel::ComputeTlimitmin()
250{
251 G4double x = (particle == positron) ?
252 0.7*msc[idx]->sqrtZ*stepmin : 0.87*msc[idx]->Z23*stepmin;
253 if(currentKinEnergy < tlow) { x *= 0.5*(1.+currentKinEnergy/tlow); }
254 return std::max(x, tlimitminfix);
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258
259#endif
260
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void flatArray(const int size, double *vect)=0
G4double GetPDGCharge() const
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double ComputeTrueStepLength(G4double geomStepLength) override
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
void StartTracking(G4Track *) override
~G4UrbanMscModel() override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeGeomPathLength(G4double truePathLength) override
G4UrbanMscModel & operator=(const G4UrbanMscModel &right)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
G4UrbanMscModel(const G4UrbanMscModel &)=delete
#define DBL_MAX
Definition: templates.hh:62