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