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