Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VMultipleScattering.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//
29// GEANT4 Class header file
30//
31//
32// File name: G4VMultipleScattering
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 12.03.2002
37//
38// Modifications:
39//
40// 16-07-03 Update GetRange interface (V.Ivanchenko)
41// 26-11-03 bugfix in AlongStepDoIt (L.Urban)
42// 25-05-04 add protection against case when range is less than steplimit (VI)
43// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
44// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
45// 15-04-05 remove boundary flag (V.Ivanchenko)
46// 07-10-05 error in a protection in GetContinuousStepLimit corrected (L.Urban)
47// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
48// 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
49// 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
50// 07-03-06 Move step limit calculation to model (V.Ivanchenko)
51// 13-05-06 Add method to access model by index (V.Ivanchenko)
52// 12-02-07 Add get/set skin (V.Ivanchenko)
53// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
54// 15-07-08 Reorder class members for further multi-thread development (VI)
55// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
56//
57// Class Description:
58//
59// It is the generic process of multiple scattering it includes common
60// part of calculations for all charged particles
61
62// -------------------------------------------------------------------
63//
64
65#ifndef G4VMultipleScattering_h
66#define G4VMultipleScattering_h 1
67
69#include "globals.hh"
70#include "G4Material.hh"
72#include "G4Track.hh"
73#include "G4Step.hh"
74#include "G4EmModelManager.hh"
75#include "G4VMscModel.hh"
76#include "G4EmParameters.hh"
77#include "G4MscStepLimitType.hh"
78
82class G4SafetyHelper;
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87{
88public:
89
90 explicit G4VMultipleScattering(const G4String& name = "msc",
92
93 ~G4VMultipleScattering() override;
94
95 //------------------------------------------------------------------------
96 // Virtual methods to be implemented for the concrete model
97 //------------------------------------------------------------------------
98
99 void ProcessDescription(std::ostream& outFile) const override;
100
101 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
102
103 // Print out of generic class parameters
104 void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&,
105 G4bool rst = false) const;
106
107protected:
108
109 virtual void StreamProcessInfo(std::ostream&) const {};
110
111public:
112
113 //------------------------------------------------------------------------
114 // Generic methods common to all ContinuousDiscrete processes
115 //------------------------------------------------------------------------
116
117 // Initialise for build of tables
118 void PreparePhysicsTable(const G4ParticleDefinition&) override;
119
120 // Build physics table during initialisation
121 void BuildPhysicsTable(const G4ParticleDefinition&) override;
122
123 // Store PhysicsTable in a file.
124 // Return false in case of failure at I/O
126 const G4String& directory,
127 G4bool ascii = false) override;
128
129 // Retrieve Physics from a file.
130 // (return true if the Physics Table can be build by using file)
131 // (return false if the process has no functionality or in case of failure)
132 // File name should is constructed as processName+particleName and the
133 // should be placed under the directory specified by the argument.
135 const G4String& directory,
136 G4bool ascii) override;
137
138 // This is called in the beginning of tracking for a new track
139 void StartTracking(G4Track*) override;
140
141 // The function overloads the corresponding function of the base
142 // class.It limits the step near to boundaries only
143 // and invokes the method GetMscContinuousStepLimit at every step.
145 const G4Track&,
146 G4double previousStepSize,
147 G4double currentMinimalStep,
148 G4double& currentSafety,
149 G4GPILSelection* selection) override;
150
151 // The function overloads the corresponding function of the base
152 // class.
154 const G4Track&,
155 G4double previousStepSize,
156 G4ForceCondition* condition) override;
157
158 // Along step actions
159 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
160
161 // This method does not used for tracking, it is intended only for tests
163 G4double previousStepSize,
164 G4double currentMinimalStep,
165 G4double& currentSafety);
166
167 // hide assignment operator
170
171 //------------------------------------------------------------------------
172 // Specific methods to set, access, modify models
173 //------------------------------------------------------------------------
174
175 // Select model in run time
176 inline G4VEmModel* SelectModel(G4double kinEnergy, size_t idx);
177
178public:
179
180 // Add model for region, smaller value of order defines which
181 // model will be selected for a given energy interval
182 void AddEmModel(G4int order, G4VMscModel*, const G4Region* region = nullptr);
183
184 // Assign a model to a process local list, to enable the list in run time
185 // the derived process should execute AddEmModel(..) for all such models
186 void SetEmModel(G4VMscModel*, G4int idx = 0);
187
188 // return a model from the local list
189 inline G4VMscModel* EmModel(size_t index = 0) const;
190
191 // Access to run time models
192 inline G4int NumberOfModels() const;
193
194 inline G4VMscModel* GetModelByIndex(G4int idx, G4bool ver = false) const;
195
196 //------------------------------------------------------------------------
197 // Get/Set parameters for simulation of multiple scattering
198 //------------------------------------------------------------------------
199
200 inline G4bool LateralDisplasmentFlag() const;
201
202 inline G4double Skin() const;
203
204 inline G4double RangeFactor() const;
205
206 inline G4double GeomFactor() const;
207
208 inline G4double PolarAngleLimit() const;
209
210 inline G4bool UseBaseMaterial() const;
211
212 inline G4MscStepLimitType StepLimitType() const;
213 inline void SetStepLimitType(G4MscStepLimitType val);
214
215 inline G4double LowestKinEnergy() const;
216 inline void SetLowestKinEnergy(G4double val);
217
218 inline const G4ParticleDefinition* FirstParticle() const;
219
220 //------------------------------------------------------------------------
221 // Run time methods
222 //------------------------------------------------------------------------
223
224protected:
225
226 // This method is not used for tracking, it returns mean free path value
227 G4double GetMeanFreePath(const G4Track& track,
228 G4double,
229 G4ForceCondition* condition) override;
230
231 // This method is not used for tracking, it returns step limit
233 G4double previousStepSize,
234 G4double currentMinimalStep,
235 G4double& currentSafety) override;
236
237private:
238
239 // ======== Parameters of the class fixed at construction =========
240
241 G4EmModelManager* modelManager;
242 G4LossTableManager* emManager;
243 G4EmParameters* theParameters;
244
245 // ======== Parameters of the class fixed at initialisation =======
246
247 G4SafetyHelper* safetyHelper = nullptr;
248 const G4ParticleDefinition* firstParticle = nullptr;
249 const G4ParticleDefinition* currParticle = nullptr;
250
251 std::vector<G4VMscModel*> mscModels;
252
253 G4double facrange = 0.04;
254 G4double lowestKinEnergy;
255
256 // ======== Cached values - may be state dependent ================
257
258protected:
259
261
262private:
263
264 G4ThreeVector fNewPosition;
265 G4ThreeVector fNewDirection;
266
267 G4VMscModel* currentModel = nullptr;
268 G4VEnergyLossProcess* fIonisation = nullptr;
269
270 G4double geomMin;
271 G4double minDisplacement2;
272 G4double physStepLimit = 0.0;
273 G4double tPathLength = 0.0;
274 G4double gPathLength = 0.0;
275
276 G4MscStepLimitType stepLimit = fUseSafety;
277 G4int numberOfModels = 0;
278
279 G4bool latDisplacement = true;
280 G4bool isIon = false;
281 G4bool fPositionChanged = false;
282 G4bool isActive = false;
283 G4bool baseMat = false;
284};
285
286// ======== Run time inline methods ================
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
289
290inline G4VEmModel*
291G4VMultipleScattering::SelectModel(G4double kinEnergy, size_t coupleIndex)
292{
293 return modelManager->SelectModel(kinEnergy, coupleIndex);
294}
295
296//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
297
299{
300 return latDisplacement;
301}
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
304
306{
307 return theParameters->MscSkin();
308}
309
310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311
313{
314 return facrange;
315}
316
317//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
318
320{
321 return theParameters->MscGeomFactor();
322}
323
324//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325
327{
328 return theParameters->MscThetaLimit();
329}
330
331//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
332
334{
335 return stepLimit;
336}
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
339
341{
342 theParameters->SetMscStepLimitType(val);
343}
344
345//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
346
348{
349 return lowestKinEnergy;
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353
355{
356 lowestKinEnergy = val;
357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
360
362{
363 return firstParticle;
364}
365
366//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367
369{
370 return (index < mscModels.size()) ? mscModels[index] : nullptr;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
374
376{
377 return numberOfModels;
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
381
382inline G4VMscModel*
384{
385 // static cast is possible inside this class
386 return static_cast<G4VMscModel*>(modelManager->GetModel(idx, ver));
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
392{
393 return baseMat;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
398#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4MscStepLimitType
@ fUseSafety
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4VEmModel * SelectModel(G4double energy, std::size_t index)
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
G4double MscThetaLimit() const
G4double MscGeomFactor() const
void SetMscStepLimitType(G4MscStepLimitType val)
G4double MscSkin() const
void AddEmModel(G4int order, G4VMscModel *, const G4Region *region=nullptr)
const G4ParticleDefinition * FirstParticle() const
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4ParticleChangeForMSC fParticleChange
G4VMscModel * GetModelByIndex(G4int idx, G4bool ver=false) const
void StartTracking(G4Track *) override
G4VMultipleScattering(G4VMultipleScattering &)=delete
void PreparePhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
void SetStepLimitType(G4MscStepLimitType val)
void ProcessDescription(std::ostream &outFile) const override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void SetLowestKinEnergy(G4double val)
void SetEmModel(G4VMscModel *, G4int idx=0)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
void StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
G4VMultipleScattering & operator=(const G4VMultipleScattering &right)=delete
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
virtual void StreamProcessInfo(std::ostream &) const
G4MscStepLimitType StepLimitType() const
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
G4VMscModel * EmModel(size_t index=0) const