Geant4_10
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: G4WentzelVIModel.hh 74726 2013-10-21 08:42:46Z gcosmo $
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"
61 #include "G4MaterialCutsCouple.hh"
63 
65 class G4LossTableManager;
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70 {
71 
72 public:
73 
74  G4WentzelVIModel(const G4String& nam = "WentzelVIUni");
75 
76  virtual ~G4WentzelVIModel();
77 
78  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
79 
80  void StartTracking(G4Track*);
81 
83  G4double KineticEnergy,
84  G4double AtomicNumber,
85  G4double AtomicWeight=0.,
86  G4double cut = DBL_MAX,
87  G4double emax= DBL_MAX);
88 
90  G4double safety);
91 
92  virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
93  G4double& currentMinimalStep);
94 
95  virtual G4double ComputeGeomPathLength(G4double truePathLength);
96 
97  virtual G4double ComputeTrueStepLength(G4double geomStepLength);
98 
99  // defines low energy limit on energy transfer to atomic electron
100  inline void SetFixedCut(G4double);
101 
102  // low energy limit on energy transfer to atomic electron
103  inline G4double GetFixedCut() const;
104 
105  // access to cross section class
107 
108 private:
109 
110  G4double ComputeXSectionPerVolume();
111 
112  inline void SetupParticle(const G4ParticleDefinition*);
113 
114  inline void DefineMaterial(const G4MaterialCutsCouple*);
115 
116  // hide assignment operator
117  G4WentzelVIModel & operator=(const G4WentzelVIModel &right);
119 
120  G4LossTableManager* theManager;
121  G4ParticleChangeForMSC* fParticleChange;
123 
124  const G4DataVector* currentCuts;
125 
126  G4double tlimitminfix;
127  G4double invsqrt12;
128  G4double fixedCut;
129 
130  // cache kinematics
131  G4double preKinEnergy;
132  G4double tPathLength;
133  G4double zPathLength;
134  G4double lambdaeff;
135  G4double currentRange;
136 
137  // data for single scattering mode
138  G4double xtsec;
139  std::vector<G4double> xsecn;
140  std::vector<G4double> prob;
141  G4int nelments;
142 
143  G4double numlimit;
144 
145  // cache material
146  G4int currentMaterialIndex;
147  const G4MaterialCutsCouple* currentCouple;
148  const G4Material* currentMaterial;
149 
150  // single scattering parameters
151  G4double cosThetaMin;
152  G4double cosThetaMax;
153  G4double cosTetMaxNuc;
154 
155  // projectile
156  const G4ParticleDefinition* particle;
157  G4double lowEnergyLimit;
158 
159  // flags
160  G4bool inside;
161  G4bool singleScatteringMode;
162 };
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166 
167 inline
168 void G4WentzelVIModel::DefineMaterial(const G4MaterialCutsCouple* cup)
169 {
170  if(cup != currentCouple) {
171  currentCouple = cup;
172  SetCurrentCouple(cup);
173  currentMaterial = cup->GetMaterial();
174  currentMaterialIndex = currentCouple->GetIndex();
175  }
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
180 inline void G4WentzelVIModel::SetupParticle(const G4ParticleDefinition* p)
181 {
182  // Initialise mass and charge
183  if(p != particle) {
184  particle = p;
185  wokvi->SetupParticle(p);
186  }
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190 
192 {
193  fixedCut = val;
194 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197 
199 {
200  return fixedCut;
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204 
206 {
207  return wokvi;
208 }
209 
210 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
211 
212 #endif
213 
G4double GetFixedCut() const
virtual ~G4WentzelVIModel()
void SetupParticle(const G4ParticleDefinition *)
G4WentzelVIModel(const G4String &nam="WentzelVIUni")
const char * p
Definition: xmltok.h:285
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX)
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void StartTracking(G4Track *)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
void SetFixedCut(G4double)
#define DBL_MAX
Definition: templates.hh:83
G4WentzelOKandVIxSection * GetWVICrossSection()
virtual G4double ComputeGeomPathLength(G4double truePathLength)
const G4Material * GetMaterial() const