Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4WentzelVIModel
35//
36// Author: V.Ivanchenko
37//
38// Creation date: 09.04.2008 from G4MuMscModel
39//
40// Modifications:
41// 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
42// compute cross sections and sample scattering angle
43//
44// Class Description:
45//
46// Implementation of the model of multiple scattering based on
47// G.Wentzel, Z. Phys. 40 (1927) 590.
48// H.W.Lewis, Phys Rev 78 (1950) 526.
49// J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
50// L.Urban, CERN-OPEN-2006-077.
51
52// -------------------------------------------------------------------
53//
54
55#ifndef G4WentzelVIModel_h
56#define G4WentzelVIModel_h 1
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60#include "G4VMscModel.hh"
63
66class G4NistManager;
67class G4Pow;
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
72{
73
74public:
75
76 G4WentzelVIModel(const G4String& nam = "WentzelVIUni");
77
78 virtual ~G4WentzelVIModel();
79
80 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
81
83
85 G4double KineticEnergy,
86 G4double AtomicNumber,
87 G4double AtomicWeight=0.,
88 G4double cut = DBL_MAX,
89 G4double emax= DBL_MAX);
90
92 G4double safety);
93
94 virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
95 G4double& currentMinimalStep);
96
97 virtual G4double ComputeGeomPathLength(G4double truePathLength);
98
99 virtual G4double ComputeTrueStepLength(G4double geomStepLength);
100
101private:
102
103 G4double ComputeXSectionPerVolume();
104
105 inline void SetupParticle(const G4ParticleDefinition*);
106
107 inline void DefineMaterial(const G4MaterialCutsCouple*);
108
109 // hide assignment operator
110 G4WentzelVIModel & operator=(const G4WentzelVIModel &right);
112
113 G4LossTableManager* theManager;
114 G4NistManager* fNistManager;
115 G4ParticleChangeForMSC* fParticleChange;
117 G4Pow* fG4pow;
118
119 const G4DataVector* currentCuts;
120
121 G4double tlimitminfix;
122 G4double invsqrt12;
123
124 // cache kinematics
125 G4double preKinEnergy;
126 G4double tPathLength;
127 G4double zPathLength;
128 G4double lambdaeff;
129 G4double currentRange;
130
131 // data for single scattering mode
132 G4double xtsec;
133 std::vector<G4double> xsecn;
134 std::vector<G4double> prob;
135 G4int nelments;
136
137 G4double numlimit;
138
139 // cache material
140 G4int currentMaterialIndex;
141 const G4MaterialCutsCouple* currentCouple;
142 const G4Material* currentMaterial;
143
144 // single scattering parameters
145 G4double cosThetaMin;
146 G4double cosThetaMax;
147 G4double cosTetMaxNuc;
148
149 // projectile
150 const G4ParticleDefinition* particle;
151 G4double lowEnergyLimit;
152
153 // flags
154 G4bool inside;
155 G4bool singleScatteringMode;
156};
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
160
161inline
162void G4WentzelVIModel::DefineMaterial(const G4MaterialCutsCouple* cup)
163{
164 if(cup != currentCouple) {
165 currentCouple = cup;
166 SetCurrentCouple(cup);
167 currentMaterial = cup->GetMaterial();
168 currentMaterialIndex = currentCouple->GetIndex();
169 }
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173
174inline void G4WentzelVIModel::SetupParticle(const G4ParticleDefinition* p)
175{
176 // Initialise mass and charge
177 if(p != particle) {
178 particle = p;
179 wokvi->SetupParticle(p);
180 }
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
185#endif
186
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const G4Material * GetMaterial() const
Definition: G4Pow.hh:54
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:370
void SetupParticle(const G4ParticleDefinition *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX)
virtual G4double ComputeGeomPathLength(G4double truePathLength)
void StartTracking(G4Track *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
virtual G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
virtual ~G4WentzelVIModel()
#define DBL_MAX
Definition: templates.hh:83