Geant4  10.00.p02
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 74657 2013-10-18 08:13:12Z 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 "G4MscStepLimitType.hh"
83 
86 class G4LossTableManager;
87 class G4SafetyHelper;
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92 {
93 public:
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) = 0;
105 
106  virtual void PrintInfo() = 0;
107 
108 protected:
109 
110  virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
111 
112 public:
113 
114  //------------------------------------------------------------------------
115  // Generic methods common to all ContinuousDiscrete processes
116  //------------------------------------------------------------------------
117 
118  // Initialise for build of tables
120 
121  // Build physics table during initialisation
123 
124  // Print out of generic class parameters
125  void PrintInfoDefinition();
126 
127  // Store PhysicsTable in a file.
128  // Return false in case of failure at I/O
130  const G4String& directory,
131  G4bool ascii = false);
132 
133  // Retrieve Physics from a file.
134  // (return true if the Physics Table can be build by using file)
135  // (return false if the process has no functionality or in case of failure)
136  // File name should is constructed as processName+particleName and the
137  // should be placed under the directory specifed by the argument.
139  const G4String& directory,
140  G4bool ascii);
141 
142  // This is called in the beginning of tracking for a new track
143  void StartTracking(G4Track*);
144 
145  // The function overloads the corresponding function of the base
146  // class.It limits the step near to boundaries only
147  // and invokes the method GetMscContinuousStepLimit at every step.
149  const G4Track&,
150  G4double previousStepSize,
151  G4double currentMinimalStep,
152  G4double& currentSafety,
153  G4GPILSelection* selection);
154 
155  // The function overloads the corresponding function of the base
156  // class.
158  const G4Track&,
159  G4double previousStepSize,
161 
162  // Along step actions
164 
165  // Post step actions
166  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
167 
168  // This method does not used for tracking, it is intended only for tests
169  G4double ContinuousStepLimit(const G4Track& track,
170  G4double previousStepSize,
171  G4double currentMinimalStep,
172  G4double& currentSafety);
173 
174  //------------------------------------------------------------------------
175  // Specific methods to set, access, modify models
176  //------------------------------------------------------------------------
177 
178  // Select model in run time
179  inline G4VEmModel* SelectModel(G4double kinEnergy, size_t idx);
180 
181 public:
182 
183  // Add model for region, smaller value of order defines which
184  // model will be selected for a given energy interval
185  void AddEmModel(G4int order, G4VEmModel*, const G4Region* region = 0);
186 
187  // Assign a model to a process
188  void SetEmModel(G4VMscModel*, G4int index = 1);
189 
190  // return the assigned model
191  G4VMscModel* EmModel(G4int index = 1) const;
192 
193  // Access to models by index
194  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
195 
196  //------------------------------------------------------------------------
197  // Get/Set parameters for simulation of multiple scattering
198  //------------------------------------------------------------------------
199 
201 
202  inline G4bool LateralDisplasmentFlag() const;
203  inline void SetLateralDisplasmentFlag(G4bool val);
204 
205  inline G4double Skin() const;
206  inline void SetSkin(G4double val);
207 
208  inline G4double RangeFactor() const;
209  inline void SetRangeFactor(G4double val);
210 
211  inline G4double GeomFactor() const;
212  inline void SetGeomFactor(G4double val);
213 
214  inline G4double PolarAngleLimit() const;
215  inline void SetPolarAngleLimit(G4double val);
216 
217  inline G4MscStepLimitType StepLimitType() const;
218  inline void SetStepLimitType(G4MscStepLimitType val);
219 
220  inline G4double LowestKinEnergy() const;
221  inline void SetLowestKinEnergy(G4double val);
222 
223  inline const G4ParticleDefinition* FirstParticle() const;
224 
225  //------------------------------------------------------------------------
226  // Run time methods
227  //------------------------------------------------------------------------
228 
229 protected:
230 
231  // This method is not used for tracking, it returns mean free path value
232  G4double GetMeanFreePath(const G4Track& track,
233  G4double,
234  G4ForceCondition* condition);
235 
236  // This method is not used for tracking, it returns step limit
238  G4double previousStepSize,
239  G4double currentMinimalStep,
240  G4double& currentSafety);
241 
242 private:
243 
244  // hide assignment operator
247 
248  // ======== Parameters of the class fixed at construction =========
249 
253 
254  // ======== Parameters of the class fixed at initialisation =======
255 
257 
258  std::vector<G4VMscModel*> mscModels;
260 
263 
265 
271 
274 
275  // ======== Cashed values - may be state dependent ================
276 
277 protected:
278 
281 
282 private:
283 
284  // cache
287 
291 
295 
297 };
298 
299 // ======== Run time inline methods ================
300 
301 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
302 
303 inline G4VEmModel*
304 G4VMultipleScattering::SelectModel(G4double kinEnergy, size_t coupleIndex)
305 {
306  return modelManager->SelectModel(kinEnergy, coupleIndex);
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
310 
312 {
313  return latDisplasment;
314 }
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
317 
319 {
320  latDisplasment = val;
321 }
322 
323 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324 
326 {
327  return skin;
328 }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
331 
333 {
334  if(val < 1.0) { skin = 0.0; }
335  else { skin = val; }
336 }
337 
338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
339 
341 {
342  return facrange;
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
346 
348 {
349  if(val > 0.0) facrange = val;
350 }
351 
352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353 
355 {
356  return facgeom;
357 }
358 
359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
360 
362 {
363  if(val > 0.0) facgeom = val;
364 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367 
369 {
370  return polarAngleLimit;
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
374 
376 {
377  if(val < 0.0) { polarAngleLimit = 0.0; }
378  else if(val > CLHEP::pi) { polarAngleLimit = CLHEP::pi; }
379  else { polarAngleLimit = val; }
380 }
381 
382 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
383 
385 {
386  return stepLimit;
387 }
388 
389 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390 
392 {
393  stepLimit = val;
394  if(val == fMinimal) { facrange = 0.2; }
395 }
396 
397 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
398 
400 {
401  return lowestKinEnergy;
402 }
403 
404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
405 
407 {
408  lowestKinEnergy = val;
409 }
410 
411 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412 
414 {
415  return firstParticle;
416 }
417 
418 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
419 
420 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
G4EmModelManager * modelManager
CLHEP::Hep3Vector G4ThreeVector
void SetLateralDisplasmentFlag(G4bool val)
G4String name
Definition: TRTMaterials.hh:40
const G4double pi
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition)
void SetEmModel(G4VMscModel *, G4int index=1)
const G4ParticleDefinition * FirstParticle() const
G4double PolarAngleLimit() const
int G4int
Definition: G4Types.hh:78
G4VEnergyLossProcess * fIonisation
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
const G4ParticleDefinition * currParticle
G4MscStepLimitType StepLimitType() const
G4ParticleChangeForMSC fParticleChange
G4double LowestKinEnergy() const
G4GPILSelection valueGPILSelectionMSC
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
bool G4bool
Definition: G4Types.hh:79
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
void SetPolarAngleLimit(G4double val)
Definition: G4Step.hh:76
G4LossTableManager * emManager
G4bool LateralDisplasmentFlag() const
virtual G4bool IsApplicable(const G4ParticleDefinition &p)=0
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
const G4ParticleDefinition * firstParticle
G4MscStepLimitType
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
G4VMultipleScattering & operator=(const G4VMultipleScattering &right)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4MscStepLimitType stepLimit
void SetGeomFactor(G4double val)
virtual void PrintInfo()=0
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &)
G4ForceCondition
void SetRangeFactor(G4double val)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
G4VMscModel * EmModel(G4int index=1) const
std::vector< G4VMscModel * > mscModels
void SetLowestKinEnergy(G4double val)
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4GPILSelection
void SetStepLimitType(G4MscStepLimitType val)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void PreparePhysicsTable(const G4ParticleDefinition &)
G4ProcessType
void SetIonisation(G4VEnergyLossProcess *)