Geant4-11
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
101protected:
102
103 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
104
105 virtual void StreamProcessInfo(std::ostream&) const {};
106
107public:
108
109 //------------------------------------------------------------------------
110 // Generic methods common to all ContinuousDiscrete processes
111 //------------------------------------------------------------------------
112
113 // Initialise for build of tables
114 void PreparePhysicsTable(const G4ParticleDefinition&) override;
115
116 // Build physics table during initialisation
117 void BuildPhysicsTable(const G4ParticleDefinition&) override;
118
119 // Store PhysicsTable in a file.
120 // Return false in case of failure at I/O
122 const G4String& directory,
123 G4bool ascii = false) override;
124
125 // Retrieve Physics from a file.
126 // (return true if the Physics Table can be build by using file)
127 // (return false if the process has no functionality or in case of failure)
128 // File name should is constructed as processName+particleName and the
129 // should be placed under the directory specified by the argument.
131 const G4String& directory,
132 G4bool ascii) override;
133
134 // This is called in the beginning of tracking for a new track
135 void StartTracking(G4Track*) override;
136
137 // The function overloads the corresponding function of the base
138 // class.It limits the step near to boundaries only
139 // and invokes the method GetMscContinuousStepLimit at every step.
141 const G4Track&,
142 G4double previousStepSize,
143 G4double currentMinimalStep,
144 G4double& currentSafety,
145 G4GPILSelection* selection) override;
146
147 // The function overloads the corresponding function of the base
148 // class.
150 const G4Track&,
151 G4double previousStepSize,
152 G4ForceCondition* condition) override;
153
154 // Along step actions
155 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
156
157 // Post step actions
158 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
159
160 // This method does not used for tracking, it is intended only for tests
162 G4double previousStepSize,
163 G4double currentMinimalStep,
164 G4double& currentSafety);
165
166 // hide assignment operator
169
170 //------------------------------------------------------------------------
171 // Specific methods to set, access, modify models
172 //------------------------------------------------------------------------
173
174 // Select model in run time
175 inline G4VEmModel* SelectModel(G4double kinEnergy, size_t idx);
176
177public:
178
179 // Add model for region, smaller value of order defines which
180 // model will be selected for a given energy interval
181 void AddEmModel(G4int order, G4VEmModel*, const G4Region* region = nullptr);
182
183 // Assign a model to a process local list, to enable the list in run time
184 // the derived process should execute AddEmModel(..) for all such models
185 void SetEmModel(G4VMscModel*, G4int idx = 0);
186
187 // return a model from the local list
188 inline G4VMscModel* EmModel(size_t index = 0) const;
189
190 // Access to run time models
191 inline G4int NumberOfModels() const;
192
193 inline G4VMscModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
194
195 //------------------------------------------------------------------------
196 // Get/Set parameters for simulation of multiple scattering
197 //------------------------------------------------------------------------
198
200
201 inline G4bool LateralDisplasmentFlag() const;
202
203 inline G4double Skin() const;
204
205 inline G4double RangeFactor() const;
206
207 inline G4double GeomFactor() const;
208
209 inline G4double PolarAngleLimit() const;
210
211 inline G4bool UseBaseMaterial() const;
212
213 inline G4MscStepLimitType StepLimitType() const;
214 inline void SetStepLimitType(G4MscStepLimitType val);
215
216 inline G4double LowestKinEnergy() const;
217 inline void SetLowestKinEnergy(G4double val);
218
219 inline const G4ParticleDefinition* FirstParticle() const;
220
221 //------------------------------------------------------------------------
222 // Run time methods
223 //------------------------------------------------------------------------
224
225protected:
226
227 // This method is not used for tracking, it returns mean free path value
228 G4double GetMeanFreePath(const G4Track& track,
229 G4double,
230 G4ForceCondition* condition) override;
231
232 // This method is not used for tracking, it returns step limit
234 G4double previousStepSize,
235 G4double currentMinimalStep,
236 G4double& currentSafety) override;
237
238private:
239
240 // Print out of generic class parameters
241 void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&,
242 G4bool rst = false) const;
243
244 // ======== Parameters of the class fixed at construction =========
245
249
250 // ======== Parameters of the class fixed at initialisation =======
251
255
256 std::vector<G4VMscModel*> mscModels;
257
260
261 // ======== Cached values - may be state dependent ================
262
263protected:
264
266
267private:
268
271
274
280
283
285 G4bool isIon = false;
289};
290
291// ======== Run time inline methods ================
292
293//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
294
295inline G4VEmModel*
296G4VMultipleScattering::SelectModel(G4double kinEnergy, size_t coupleIndex)
297{
298 return modelManager->SelectModel(kinEnergy, coupleIndex);
299}
300
301//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
302
304{
305 return latDisplacement;
306}
307
308//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
309
311{
312 return theParameters->MscSkin();
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316
318{
319 return facrange;
320}
321
322//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
323
325{
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
330
332{
334}
335
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
337
339{
341}
342
343//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
344
346{
348}
349
350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
351
353{
354 return lowestKinEnergy;
355}
356
357//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
358
360{
361 lowestKinEnergy = val;
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
365
367{
368 return firstParticle;
369}
370
371//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
372
374{
375 return (index < mscModels.size()) ? mscModels[index] : nullptr;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
379
381{
382 return numberOfModels;
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386
387inline G4VMscModel*
389{
390 return static_cast<G4VMscModel*>(modelManager->GetModel(idx, ver));
391}
392
393//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
394
396{
397 return baseMat;
398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401
402#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 * GetModel(G4int idx, G4bool ver=false)
G4VEmModel * SelectModel(G4double energy, size_t index)
G4double MscThetaLimit() const
G4MscStepLimitType MscStepLimitType() const
G4double MscGeomFactor() const
void SetMscStepLimitType(G4MscStepLimitType val)
G4double MscSkin() const
Definition: G4Step.hh:62
void SetIonisation(G4VEnergyLossProcess *)
G4EmModelManager * modelManager
const G4ParticleDefinition * FirstParticle() const
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4ParticleChangeForMSC fParticleChange
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
G4LossTableManager * emManager
G4VEnergyLossProcess * fIonisation
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)
G4VMscModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double PolarAngleLimit() const
void SetStepLimitType(G4MscStepLimitType val)
void ProcessDescription(std::ostream &outFile) const override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
const G4ParticleDefinition * firstParticle
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
G4MscStepLimitType stepLimit
std::vector< G4VMscModel * > mscModels
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4bool LateralDisplasmentFlag() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
virtual void StreamProcessInfo(std::ostream &) const
G4MscStepLimitType StepLimitType() const
G4double LowestKinEnergy() const
const G4ParticleDefinition * currParticle
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
G4VMscModel * EmModel(size_t index=0) const
const char * name(G4int ptype)