Geant4-11
G4VMscModel.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// GEANT4 Class header file
29//
30//
31// File name: G4VMscModel
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 07.03.2008
36//
37// Modifications:
38// 07.04.2009 V.Ivanchenko moved msc methods from G4VEmModel to G4VMscModel
39// 26.03.2012 V.Ivanchenko added transport x-section pointer and Get?Set methods
40//
41// Class Description:
42//
43// General interface to msc models
44
45// -------------------------------------------------------------------
46//
47#ifndef G4VMscModel_h
48#define G4VMscModel_h 1
49
51
52#include "G4VEmModel.hh"
53#include "G4MscStepLimitType.hh"
54#include "globals.hh"
55#include "G4ThreeVector.hh"
56#include "G4Track.hh"
57#include "G4SafetyHelper.hh"
58#include "G4PhysicsTable.hh"
59#include "G4ThreeVector.hh"
60#include <vector>
61
65
66class G4VMscModel : public G4VEmModel
67{
68
69public:
70
71 explicit G4VMscModel(const G4String& nam);
72
73 ~G4VMscModel() override;
74
76 G4double& stepLimit) = 0;
77
78 virtual G4double ComputeGeomPathLength(G4double truePathLength) = 0;
79
80 virtual G4double ComputeTrueStepLength(G4double geomPathLength) = 0;
81
83 G4double safety) = 0;
84
86
87 void DumpParameters(std::ostream& out) const;
88
89 // empty method
90 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
92 const G4DynamicParticle*,
93 G4double tmin, G4double tmax) override;
94
95 //================================================================
96 // Set parameters of multiple scattering models
97 //================================================================
98
100
101 inline void SetLateralDisplasmentFlag(G4bool val);
102
103 inline void SetRangeFactor(G4double);
104
105 inline void SetGeomFactor(G4double);
106
107 inline void SetSkin(G4double);
108
109 inline void SetLambdaLimit(G4double);
110
111 inline void SetSafetyFactor(G4double);
112
113 inline void SetSampleZ(G4bool);
114
115 //================================================================
116 // Get/Set access to Physics Tables
117 //================================================================
118
119 inline G4VEnergyLossProcess* GetIonisation() const;
120
122 const G4ParticleDefinition* part);
123
124 //================================================================
125 // Run time methods
126 //================================================================
127
128protected:
129
130 // initialisation of the ParticleChange for the model
131 // initialisation of interface with geometry and ionisation
134
135 // convert true length to geometry length
136 inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
137
138public:
139
140 // compute safety
142 G4double limit= DBL_MAX);
143
144 // compute linear distance to a geometry boundary
145 inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
146 G4double limit);
147
149 G4double kineticEnergy,
150 const G4MaterialCutsCouple* couple);
151
153 G4double kineticEnergy,
154 const G4MaterialCutsCouple* couple,
155 G4double logKineticEnergy);
156
158 G4double kineticEnergy,
159 const G4MaterialCutsCouple* couple);
160
162 G4double kineticEnergy,
163 const G4MaterialCutsCouple* couple,
164 G4double logKineticEnergy);
165
167 G4double range,
168 const G4MaterialCutsCouple* couple);
169
170 // G4MaterialCutsCouple should be defined before call to this method
171 inline
173 G4double kinEnergy);
174
175 inline
177 G4double kinEnergy,
178 G4double logKinEnergy);
179
180 // hide assignment operator
181 G4VMscModel & operator=(const G4VMscModel &right) = delete;
182 G4VMscModel(const G4VMscModel&) = delete;
183
184private:
185
189
193
194protected:
195
204
207
210
211};
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217{
218 if(!IsLocked()) { latDisplasment = val; }
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222
224{
225 if(!IsLocked()) { skin = val; }
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
229
231{
232 if(!IsLocked()) { facrange = val; }
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236
238{
239 if(!IsLocked()) { facgeom = val; }
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243
245{
246 if(!IsLocked()) { lambdalimit = val; }
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
252{
253 if(!IsLocked()) { facsafety = val; }
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257
259{
260 if(!IsLocked()) { steppingAlgorithm = val; }
261}
262
263//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264
266{
267 if(!IsLocked()) { samplez = val; }
268}
269
270//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271
273 G4double limit)
274{
275 return safetyHelper->ComputeSafety(position, limit);
276}
277
278//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
279
281 G4double& glength)
282{
283 glength = ComputeGeomPathLength(tlength);
284 // should return true length
285 return tlength;
286}
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289
291 G4double& presafety,
292 G4double limit)
293{
295 track.GetStep()->GetPreStepPoint()->GetPosition(),
296 track.GetMomentumDirection(),
297 limit, presafety);
298}
299
300//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
301
303{
304 return ionisation;
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
308
310 const G4ParticleDefinition* part)
311{
312 ionisation = p;
313 currentPart = part;
314}
315
316//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
317
318inline G4double
320 G4double ekin)
321{
322 G4double x;
323 if (nullptr != xSectionTable) {
324 x = pFactor*(*xSectionTable)[basedCoupleIndex]->Value(ekin)/(ekin*ekin);
325 } else {
326 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
327 }
328 return (x > 0.0) ? 1.0/x : DBL_MAX;
329}
330
331inline G4double
333 G4double ekin, G4double logekin)
334{
335 G4double x;
336 if (nullptr != xSectionTable) {
337 x = pFactor*(*xSectionTable)[basedCoupleIndex]->LogVectorValue(ekin, logekin)/(ekin*ekin);
338 } else {
339 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
340 }
341 return (x > 0.0) ? 1.0/x : DBL_MAX;
342}
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345
346#endif
G4MscStepLimitType
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
const G4Step * GetStep() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:426
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:427
G4double pFactor
Definition: G4VEmModel.hh:432
G4bool IsLocked() const
Definition: G4VEmModel.hh:877
size_t basedCoupleIndex
Definition: G4VEmModel.hh:449
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
G4double dtrl
Definition: G4VMscModel.hh:200
void SetLambdaLimit(G4double)
Definition: G4VMscModel.hh:244
void DumpParameters(std::ostream &out) const
Definition: G4VMscModel.cc:137
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:159
void SetSkin(G4double)
Definition: G4VMscModel.hh:223
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:230
G4double facrange
Definition: G4VMscModel.hh:196
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:290
G4bool samplez
Definition: G4VMscModel.hh:205
G4double skin
Definition: G4VMscModel.hh:199
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:309
G4double localtkin
Definition: G4VMscModel.hh:191
G4double geomMin
Definition: G4VMscModel.hh:202
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:319
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:216
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:59
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:78
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:187
G4double dedx
Definition: G4VMscModel.hh:190
void SetSafetyFactor(G4double)
Definition: G4VMscModel.hh:251
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double localrange
Definition: G4VMscModel.hh:192
G4VEnergyLossProcess * GetIonisation() const
Definition: G4VMscModel.hh:302
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:224
G4double geomMax
Definition: G4VMscModel.hh:203
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:189
G4double lambdalimit
Definition: G4VMscModel.hh:201
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:209
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:186
G4VMscModel(const G4VMscModel &)=delete
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:237
~G4VMscModel() override
Definition: G4VMscModel.cc:72
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
Definition: G4VMscModel.cc:151
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:280
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4VMscModel & operator=(const G4VMscModel &right)=delete
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:272
G4double facsafety
Definition: G4VMscModel.hh:198
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:208
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:116
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:258
G4double facgeom
Definition: G4VMscModel.hh:197
virtual G4double ComputeGeomPathLength(G4double truePathLength)=0
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
const G4ParticleDefinition * currentPart
Definition: G4VMscModel.hh:188
void SetSampleZ(G4bool)
Definition: G4VMscModel.hh:265
#define DBL_MAX
Definition: templates.hh:62