Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelVIModel.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: G4WentzelVIModel
34//
35// Author: V.Ivanchenko
36//
37// Creation date: 09.04.2008 from G4MuMscModel
38//
39// Modifications:
40// 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
41// compute cross sections and sample scattering angle
42//
43// Class Description:
44//
45// Implementation of the model of multiple scattering based on
46// G.Wentzel, Z. Phys. 40 (1927) 590.
47// H.W.Lewis, Phys Rev 78 (1950) 526.
48// J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
49// L.Urban, CERN-OPEN-2006-077.
50
51// -------------------------------------------------------------------
52//
53
54#ifndef G4WentzelVIModel_h
55#define G4WentzelVIModel_h 1
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
59#include "G4VMscModel.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
66{
67
68public:
69
70 explicit G4WentzelVIModel(G4bool comb=true, const G4String& nam = "WentzelVIUni");
71
72 ~G4WentzelVIModel() override;
73
74 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
75
77 G4VEmModel* masterModel) override;
78
79 void StartTracking(G4Track*) override;
80
82 G4double KineticEnergy,
83 G4double AtomicNumber,
84 G4double AtomicWeight=0.,
85 G4double cut = DBL_MAX,
86 G4double emax= DBL_MAX) override;
87
89 G4double safety) override;
90
93 G4double& currentMinimalStep) override;
94
95 G4double ComputeGeomPathLength(G4double truePathLength) override;
96
97 G4double ComputeTrueStepLength(G4double geomStepLength) override;
98
99 // defines low energy limit on energy transfer to atomic electron
100 void SetFixedCut(G4double);
101
102 // low energy limit on energy transfer to atomic electron
103 G4double GetFixedCut() const;
104
105 // access to cross section class
107
109
111
112 G4bool UseSecondMoment() const;
113
115
118 G4double kineticEnergy);
119
121
123
126
127protected:
128
130
131 inline void SetupParticle(const G4ParticleDefinition*);
132
133private:
134
135 G4double ComputeSecondMoment(const G4ParticleDefinition*,
136 G4double kineticEnergy);
137
138protected:
139
142 const G4Material* currentMaterial = nullptr;
143
146 const G4DataVector* currentCuts = nullptr;
148
153
154 // cache kinematics
161
163
164 // cache kinematics
166
167 // single scattering parameters
171
173 size_t idx2 = 0;
174
175 // data for single scattering mode
177
178 // flags
182
183 std::vector<G4double> xsecn;
184 std::vector<G4double> prob;
185};
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189
191{
192 // Initialise mass and charge
193 if(p != particle) {
194 particle = p;
196 }
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200
202{
203 fixedCut = val;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207
209{
210 return fixedCut;
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214
216{
217 if(ptr != wokvi) {
218 delete wokvi;
219 wokvi = ptr;
220 }
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
226{
227 return wokvi;
228}
229
230//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231
233{
234 useSecondMoment = val;
235}
236
237//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238
240{
241 return useSecondMoment;
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245
247{
248 return fSecondMoments;
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252
253inline G4double
255 const G4MaterialCutsCouple* couple,
256 G4double ekin)
257{
258 G4double x = 0.0;
259 if(useSecondMoment) {
260 DefineMaterial(couple);
261 x = (fSecondMoments) ?
262 (*fSecondMoments)[(*theDensityIdx)[currentMaterialIndex]]->Value(ekin, idx2)
263 *(*theDensityFactor)[currentMaterialIndex]/(ekin*ekin)
264 : ComputeSecondMoment(part, ekin);
265 }
266 return x;
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270
271#endif
272
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetupParticle(const G4ParticleDefinition *)
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4WentzelVIModel & operator=(const G4WentzelVIModel &right)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
G4double ComputeTransportXSectionPerVolume(G4double cosTheta)
G4bool UseSecondMoment() const
const G4MaterialCutsCouple * currentCouple
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
G4PhysicsTable * GetSecondMomentTable()
std::vector< G4double > prob
void SetFixedCut(G4double)
G4WentzelVIModel(const G4WentzelVIModel &)=delete
G4double SecondMoment(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
~G4WentzelVIModel() override
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ComputeTrueStepLength(G4double geomStepLength) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double GetFixedCut() const
G4PhysicsTable * fSecondMoments
void SetWVICrossSection(G4WentzelOKandVIxSection *)
G4WentzelOKandVIxSection * GetWVICrossSection()
std::vector< G4double > xsecn
G4double ComputeGeomPathLength(G4double truePathLength) override
G4ParticleChangeForMSC * fParticleChange
const G4ParticleDefinition * particle
void SetUseSecondMoment(G4bool)
const G4DataVector * currentCuts
const G4Material * currentMaterial
void SetupParticle(const G4ParticleDefinition *)
void StartTracking(G4Track *) override
void SetSingleScatteringFactor(G4double)
G4WentzelOKandVIxSection * wokvi
#define DBL_MAX
Definition: templates.hh:62