Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4WentzelVIRelModel.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: G4WentzelVIRelModel.hh 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4WentzelVIRelModel
35 //
36 // Author: V.Ivanchenko
37 //
38 // Creation date: 08.06.2012 from G4WentzelVIModel
39 //
40 // Modifications:
41 //
42 // Class Description:
43 //
44 // Implementation of the model of multiple scattering based on
45 // G.Wentzel, Z. Phys. 40 (1927) 590.
46 // H.W.Lewis, Phys Rev 78 (1950) 526.
47 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
48 // L.Urban, CERN-OPEN-2006-077.
49 
50 // -------------------------------------------------------------------
51 //
52 
53 #ifndef G4WentzelVIRelModel_h
54 #define G4WentzelVIRelModel_h 1
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
58 #include "G4VMscModel.hh"
59 #include "G4MaterialCutsCouple.hh"
61 
63 class G4LossTableManager;
64 class G4NistManager;
65 class G4Pow;
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70 {
71 
72 public:
73 
74  explicit G4WentzelVIRelModel(G4bool combined = true);
75 
76  virtual ~G4WentzelVIRelModel();
77 
78  virtual void Initialise(const G4ParticleDefinition*,
79  const G4DataVector&) override;
80 
81  virtual void StartTracking(G4Track*) override;
82 
84  G4double KineticEnergy,
85  G4double AtomicNumber,
86  G4double AtomicWeight=0.,
87  G4double cut = DBL_MAX,
88  G4double emax= DBL_MAX) override;
89 
91  G4double safety) override;
92 
93  virtual G4double
95  G4double& currentMinimalStep) override;
96 
97  virtual G4double ComputeGeomPathLength(G4double truePathLength) override;
98 
99  virtual G4double ComputeTrueStepLength(G4double geomStepLength) override;
100 
101 private:
102 
103  G4double ComputeXSectionPerVolume();
104 
105  inline void SetupParticle(const G4ParticleDefinition*);
106 
107  inline void DefineMaterial(const G4MaterialCutsCouple*);
108 
109  // hide assignment operator
110  G4WentzelVIRelModel & operator=(const G4WentzelVIRelModel &right) = delete;
111  G4WentzelVIRelModel(const G4WentzelVIRelModel&) = delete;
112 
113  G4LossTableManager* theManager;
114  G4NistManager* fNistManager;
115  G4ParticleChangeForMSC* fParticleChange;
116  G4WentzelVIRelXSection* wokvi;
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 isCombined;
155  G4bool inside;
156  G4bool singleScatteringMode;
157 };
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161 
162 inline
163 void G4WentzelVIRelModel::DefineMaterial(const G4MaterialCutsCouple* cup)
164 {
165  if(cup != currentCouple) {
166  currentCouple = cup;
167  SetCurrentCouple(cup);
168  currentMaterial = cup->GetMaterial();
169  currentMaterialIndex = currentCouple->GetIndex();
170  }
171 }
172 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
175 inline void G4WentzelVIRelModel::SetupParticle(const G4ParticleDefinition* p)
176 {
177  // Initialise mass and charge
178  if(p != particle) {
179  particle = p;
180  wokvi->SetupParticle(p);
181  }
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
186 #endif
187 
Definition: G4Pow.hh:56
const char * p
Definition: xmltok.h:285
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
virtual void StartTracking(G4Track *) override
int G4int
Definition: G4Types.hh:78
void SetupParticle(const G4ParticleDefinition *)
bool G4bool
Definition: G4Types.hh:79
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
G4WentzelVIRelModel(G4bool combined=true)
static const G4double emax
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
const G4Material * GetMaterial() const