Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
138 // should be set before initialisation
139 inline void SetUseSplineForMSC(G4bool val);
140
141public:
142
143 // compute safety
145 G4double limit= DBL_MAX);
146
147 // compute linear distance to a geometry boundary
148 inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
149 G4double limit);
150
152 G4double kineticEnergy,
153 const G4MaterialCutsCouple* couple);
154
156 G4double kineticEnergy,
157 const G4MaterialCutsCouple* couple,
158 G4double logKineticEnergy);
159
161 G4double kineticEnergy,
162 const G4MaterialCutsCouple* couple);
163
165 G4double kineticEnergy,
166 const G4MaterialCutsCouple* couple,
167 G4double logKineticEnergy);
168
170 G4double range,
171 const G4MaterialCutsCouple* couple);
172
173 // G4MaterialCutsCouple should be defined before call to this method
174 inline
176 G4double kinEnergy);
177
178 inline
180 G4double kinEnergy,
181 G4double logKinEnergy);
182
183 // hide assignment operator
184 G4VMscModel & operator=(const G4VMscModel &right) = delete;
185 G4VMscModel(const G4VMscModel&) = delete;
186
187private:
188
189 G4SafetyHelper* safetyHelper = nullptr;
190 G4VEnergyLossProcess* ionisation = nullptr;
191 const G4ParticleDefinition* currentPart = nullptr;
192
193 G4double dedx = 0.0;
194 G4double localtkin = 0.0;
195 G4double localrange = DBL_MAX;
196
197protected:
198
207
210
213
214private:
215
216 G4bool useSpline = true;
217};
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
221
223{
224 if(!IsLocked()) { latDisplasment = val; }
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
230{
231 if(!IsLocked()) { skin = val; }
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235
237{
238 if(!IsLocked()) { facrange = val; }
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242
244{
245 if(!IsLocked()) { facgeom = val; }
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
249
251{
252 if(!IsLocked()) { lambdalimit = val; }
253}
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
256
258{
259 if(!IsLocked()) { facsafety = val; }
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263
265{
266 if(!IsLocked()) { steppingAlgorithm = val; }
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270
272{
273 if(!IsLocked()) { samplez = val; }
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
277
279 G4double limit)
280{
281 return safetyHelper->ComputeSafety(position, limit);
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285
287 G4double& glength)
288{
289 glength = ComputeGeomPathLength(tlength);
290 // should return true length
291 return tlength;
292}
293
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
295
297 G4double& presafety,
298 G4double limit)
299{
300 return safetyHelper->CheckNextStep(
301 track.GetStep()->GetPreStepPoint()->GetPosition(),
302 track.GetMomentumDirection(),
303 limit, presafety);
304}
305
306//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
307
309{
310 return ionisation;
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
314
316 const G4ParticleDefinition* part)
317{
318 ionisation = p;
319 currentPart = part;
320}
321
322//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
323
324inline G4double
326 G4double ekin)
327{
328 G4double x;
329 if (nullptr != xSectionTable) {
330 x = pFactor*(*xSectionTable)[basedCoupleIndex]->Value(ekin)/(ekin*ekin);
331 } else {
332 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
333 }
334 return (x > 0.0) ? 1.0/x : DBL_MAX;
335}
336
337//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
338
339inline G4double
341 G4double ekin, G4double logekin)
342{
343 G4double x;
344 if (nullptr != xSectionTable) {
345 x = pFactor*(*xSectionTable)[basedCoupleIndex]->LogVectorValue(ekin, logekin)/(ekin*ekin);
346 } else {
347 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
348 }
349 return (x > 0.0) ? 1.0/x : DBL_MAX;
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
353
355{
356 useSpline = val;
357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
360
361#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:422
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:181
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:423
G4double pFactor
Definition: G4VEmModel.hh:428
G4bool IsLocked() const
Definition: G4VEmModel.hh:856
size_t basedCoupleIndex
Definition: G4VEmModel.hh:445
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
G4double dtrl
Definition: G4VMscModel.hh:203
void SetLambdaLimit(G4double)
Definition: G4VMscModel.hh:250
void DumpParameters(std::ostream &out) const
Definition: G4VMscModel.cc:136
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:158
void SetSkin(G4double)
Definition: G4VMscModel.hh:229
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:236
G4double facrange
Definition: G4VMscModel.hh:199
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:296
G4bool samplez
Definition: G4VMscModel.hh:211
G4double skin
Definition: G4VMscModel.hh:202
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:315
G4double geomMin
Definition: G4VMscModel.hh:205
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:325
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:222
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:77
void SetSafetyFactor(G4double)
Definition: G4VMscModel.hh:257
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4VEnergyLossProcess * GetIonisation() const
Definition: G4VMscModel.hh:308
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:223
G4double geomMax
Definition: G4VMscModel.hh:206
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:188
G4double lambdalimit
Definition: G4VMscModel.hh:204
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:209
G4VMscModel(const G4VMscModel &)=delete
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:243
~G4VMscModel() override
void SetUseSplineForMSC(G4bool val)
Definition: G4VMscModel.hh:354
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
Definition: G4VMscModel.cc:150
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:286
G4bool latDisplasment
Definition: G4VMscModel.hh:212
G4VMscModel & operator=(const G4VMscModel &right)=delete
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:278
G4double facsafety
Definition: G4VMscModel.hh:201
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:208
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:115
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:264
G4double facgeom
Definition: G4VMscModel.hh:200
virtual G4double ComputeGeomPathLength(G4double truePathLength)=0
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
void SetSampleZ(G4bool)
Definition: G4VMscModel.hh:271
#define DBL_MAX
Definition: templates.hh:62