source: trunk/source/processes/electromagnetic/utils/include/G4VMultipleScattering.hh @ 1340

Last change on this file since 1340 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 16.5 KB
Line 
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,v 1.63 2010/03/10 18:29:51 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4VMultipleScattering
35//
36// Author:        Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 12.03.2002
39//
40// Modifications:
41//
42// 16-07-03 Update GetRange interface (V.Ivanchenko)
43//
44//
45// Class Description:
46//
47// It is the generic process of multiple scattering it includes common
48// part of calculations for all charged particles
49//
50// 26-11-03 bugfix in AlongStepDoIt (L.Urban)
51// 25-05-04 add protection against case when range is less than steplimit (VI)
52// 30-06-04 make destructor virtual (V.Ivanchenko)
53// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
54// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
55// 15-04-05 optimize internal interfaces (V.Ivanchenko)
56// 15-04-05 remove boundary flag (V.Ivanchenko)
57// 07-10-05 error in a protection in GetContinuousStepLimit corrected (L.Urban)
58// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
59// 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
60// 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
61// 07-03-06 Move step limit calculation to model (V.Ivanchenko)
62// 13-05-06 Add method to access model by index (V.Ivanchenko)
63// 12-02-07 Add get/set skin (V.Ivanchenko)
64// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
65// 15-07-08 Reorder class members for further multi-thread development (VI)
66// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
67//
68
69// -------------------------------------------------------------------
70//
71
72#ifndef G4VMultipleScattering_h
73#define G4VMultipleScattering_h 1
74
75#include "G4VContinuousDiscreteProcess.hh"
76#include "globals.hh"
77#include "G4Material.hh"
78#include "G4MaterialCutsCouple.hh"
79#include "G4ParticleChangeForMSC.hh"
80#include "G4Track.hh"
81#include "G4Step.hh"
82#include "G4EmModelManager.hh"
83#include "G4VMscModel.hh"
84#include "G4MscStepLimitType.hh"
85
86class G4ParticleDefinition;
87class G4DataVector;
88class G4PhysicsTable;
89class G4PhysicsVector;
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93class G4VMultipleScattering : public G4VContinuousDiscreteProcess
94{
95public:
96
97  G4VMultipleScattering(const G4String& name = "msc",
98                        G4ProcessType type = fElectromagnetic);
99
100  virtual ~G4VMultipleScattering();
101
102  //------------------------------------------------------------------------
103  // Virtual methods to be implemented for the concrete model
104  //------------------------------------------------------------------------
105
106  virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
107
108  virtual void PrintInfo() = 0;
109
110protected:
111
112  virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
113
114public:
115
116  //------------------------------------------------------------------------
117  // Generic methods common to all ContinuousDiscrete processes
118  //------------------------------------------------------------------------
119
120  // Initialise for build of tables
121  void PreparePhysicsTable(const G4ParticleDefinition&);
122 
123  // Build physics table during initialisation
124  void BuildPhysicsTable(const G4ParticleDefinition&);
125
126  // Print out of generic class parameters
127  void PrintInfoDefinition();
128
129  // Store PhysicsTable in a file.
130  // Return false in case of failure at I/O
131  G4bool StorePhysicsTable(const G4ParticleDefinition*,
132                           const G4String& directory,
133                           G4bool ascii = false);
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 specifed by the argument.
140  G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
141                              const G4String& directory,
142                              G4bool ascii);
143
144  // The function overloads the corresponding function of the base
145  // class.It limits the step near to boundaries only
146  // and invokes the method GetMscContinuousStepLimit at every step.
147  G4double AlongStepGetPhysicalInteractionLength(
148                                            const G4Track&,
149                                            G4double  previousStepSize,
150                                            G4double  currentMinimalStep,
151                                            G4double& currentSafety,
152                                            G4GPILSelection* selection);
153
154  // The function overloads the corresponding function of the base
155  // class.
156  G4double PostStepGetPhysicalInteractionLength(
157                                            const G4Track&,
158                                            G4double  previousStepSize,
159                                            G4ForceCondition* condition);
160
161  // Along step actions
162  G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
163
164  // Post step actions
165  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
166
167  // This method does not used for tracking, it is intended only for tests
168  G4double ContinuousStepLimit(const G4Track& track,
169                               G4double previousStepSize,
170                               G4double currentMinimalStep,
171                               G4double& currentSafety);
172
173  //------------------------------------------------------------------------
174  // Specific methods to build and access Physics Tables
175  //------------------------------------------------------------------------
176
177  // Build empty Physics Vector
178  G4PhysicsVector* PhysicsVector(const G4MaterialCutsCouple*);
179
180  inline void SetBinning(G4int nbins);
181  inline G4int Binning() const;
182
183  inline void SetMinKinEnergy(G4double e);
184  inline G4double MinKinEnergy() const;
185
186  inline void SetMaxKinEnergy(G4double e);
187  inline G4double MaxKinEnergy() const;
188
189  inline void SetBuildLambdaTable(G4bool val);
190
191  inline G4PhysicsTable* LambdaTable() const;
192
193  // access particle type
194  inline const G4ParticleDefinition* Particle() const;
195
196  //------------------------------------------------------------------------
197  // Specific methods to set, access, modify models
198  //------------------------------------------------------------------------
199
200protected:
201  // Select model in run time
202  inline G4VEmModel* SelectModel(G4double kinEnergy);
203
204public:
205  // Select model in run time
206  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 
207                                            size_t& idxRegion) const;
208
209  // Add model for region, smaller value of order defines which
210  // model will be selected for a given energy interval 
211  void AddEmModel(G4int order, G4VEmModel*, const G4Region* region = 0);
212
213  // Assign a model to a process
214  void SetModel(G4VMscModel*, G4int index = 1);
215 
216  // return the assigned model
217  G4VMscModel* Model(G4int index = 1);
218
219  // Access to models by index
220  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
221
222  //------------------------------------------------------------------------
223  // Get/Set parameters for simulation of multiple scattering
224  //------------------------------------------------------------------------
225
226  inline G4bool LateralDisplasmentFlag() const;
227  inline void SetLateralDisplasmentFlag(G4bool val);
228
229  inline G4double Skin() const;
230  inline void SetSkin(G4double val);
231
232  inline G4double RangeFactor() const;
233  inline void SetRangeFactor(G4double val);
234
235  inline G4double GeomFactor() const;
236  inline void SetGeomFactor(G4double val);
237
238  inline G4double PolarAngleLimit() const;
239  inline void SetPolarAngleLimit(G4double val);
240
241  inline G4MscStepLimitType StepLimitType() const;
242  inline void SetStepLimitType(G4MscStepLimitType val);
243
244  //------------------------------------------------------------------------
245  // Run time methods
246  //------------------------------------------------------------------------
247
248protected:
249
250  // This method is not used for tracking, it returns mean free path value
251  G4double GetMeanFreePath(const G4Track& track,
252                           G4double,
253                           G4ForceCondition* condition);
254
255  // This method is not used for tracking, it returns step limit
256  G4double GetContinuousStepLimit(const G4Track& track,
257                                  G4double previousStepSize,
258                                  G4double currentMinimalStep,
259                                  G4double& currentSafety);
260
261  // This method returns inversed transport cross section
262  inline G4double GetLambda(const G4ParticleDefinition* p, 
263                            G4double& kineticEnergy);
264
265  // defines current material in run time
266  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
267
268  inline const G4MaterialCutsCouple* CurrentMaterialCutsCouple() const; 
269
270private:
271
272  // hide  assignment operator
273  G4VMultipleScattering(G4VMultipleScattering &);
274  G4VMultipleScattering & operator=(const G4VMultipleScattering &right);
275
276  // ======== Parameters of the class fixed at construction =========
277
278  G4EmModelManager*           modelManager;
279  G4bool                      buildLambdaTable;
280
281  // ======== Parameters of the class fixed at initialisation =======
282
283  std::vector<G4VMscModel*>   mscModels;
284
285  G4PhysicsTable*             theLambdaTable;
286  const G4ParticleDefinition* firstParticle;
287
288  G4MscStepLimitType          stepLimit;
289
290  G4double                    minKinEnergy;
291  G4double                    maxKinEnergy;
292  G4double                    skin;
293  G4double                    facrange;
294  G4double                    facgeom;
295  G4double                    polarAngleLimit;
296
297  G4int                       nBins;
298
299  G4bool                      latDisplasment;
300  G4bool                      isIon;
301
302  // ======== Cashed values - may be state dependent ================
303
304protected:
305
306  G4GPILSelection             valueGPILSelectionMSC;
307  G4ParticleChangeForMSC      fParticleChange;
308
309private:
310
311  G4VMscModel*                currentModel;
312
313  // cache
314  const G4ParticleDefinition* currentParticle;
315  const G4MaterialCutsCouple* currentCouple;
316  size_t                      currentMaterialIndex;
317
318};
319
320// ======== Run time inline methods ================
321
322inline const G4MaterialCutsCouple* 
323G4VMultipleScattering::CurrentMaterialCutsCouple() const
324{
325  return currentCouple;
326} 
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329
330inline 
331void G4VMultipleScattering::DefineMaterial(const G4MaterialCutsCouple* couple)
332{
333  if(couple != currentCouple) {
334    currentCouple   = couple;
335    currentMaterialIndex = couple->GetIndex();
336  }
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340
341inline 
342G4double G4VMultipleScattering::GetLambda(const G4ParticleDefinition* p, 
343                                          G4double& e)
344{
345  G4double x;
346  if(theLambdaTable) {
347    x = ((*theLambdaTable)[currentMaterialIndex])->Value(e);
348  } else {
349    x = currentModel->CrossSection(currentCouple,p,e);
350  }
351  if(x > DBL_MIN) { x = 1./x; }
352  else            { x = DBL_MAX; } 
353  return x;
354}
355
356//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
357
358inline G4VEmModel* G4VMultipleScattering::SelectModel(G4double kinEnergy)
359{
360  return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
361}
362
363//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
364
365inline G4VEmModel* G4VMultipleScattering::SelectModelForMaterial(
366                   G4double kinEnergy, size_t& idxRegion) const
367{
368  return modelManager->SelectModel(kinEnergy, idxRegion);
369}
370
371// ======== Get/Set inline methods used at initialisation ================
372
373inline void G4VMultipleScattering::SetBinning(G4int nbins)
374{
375  nBins = nbins;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
379
380inline G4int G4VMultipleScattering::Binning() const
381{
382  return nBins;
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386
387inline void G4VMultipleScattering::SetMinKinEnergy(G4double e)
388{
389  minKinEnergy = e;
390}
391
392//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
393
394inline G4double G4VMultipleScattering::MinKinEnergy() const
395{
396  return minKinEnergy;
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
400
401inline void G4VMultipleScattering::SetMaxKinEnergy(G4double e)
402{
403  maxKinEnergy = e;
404}
405
406//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
407
408inline G4double G4VMultipleScattering::MaxKinEnergy() const
409{
410  return maxKinEnergy;
411}
412
413//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
414
415inline  void G4VMultipleScattering::SetBuildLambdaTable(G4bool val)
416{
417  buildLambdaTable = val;
418}
419
420//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
421
422inline G4PhysicsTable* G4VMultipleScattering::LambdaTable() const
423{
424  return theLambdaTable;
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
429inline  const G4ParticleDefinition* G4VMultipleScattering::Particle() const
430{
431  return currentParticle;
432}
433
434//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
435
436inline  G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
437{
438  return latDisplasment;
439}
440
441//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
442
443inline  void G4VMultipleScattering::SetLateralDisplasmentFlag(G4bool val)
444{
445  latDisplasment = val;
446}
447
448//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
449
450inline  G4double G4VMultipleScattering::Skin() const
451{
452  return skin;
453}
454
455//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
456
457inline  void G4VMultipleScattering::SetSkin(G4double val)
458{
459  if(val < 1.0) { skin = 0.0; }
460  else          { skin = val; }
461}
462
463//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
464
465inline  G4double G4VMultipleScattering::RangeFactor() const
466{
467  return facrange;
468}
469
470//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
471
472inline  void G4VMultipleScattering::SetRangeFactor(G4double val)
473{
474  if(val > 0.0) facrange = val;
475}
476
477//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
478
479inline  G4double G4VMultipleScattering::GeomFactor() const
480{
481  return facgeom;
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
486inline  void G4VMultipleScattering::SetGeomFactor(G4double val)
487{
488  if(val > 0.0) facgeom = val;
489}
490
491//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492
493inline  G4double G4VMultipleScattering::PolarAngleLimit() const
494{
495  return polarAngleLimit;
496}
497
498//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
499
500inline  void G4VMultipleScattering::SetPolarAngleLimit(G4double val)
501{
502  if(val < 0.0)            { polarAngleLimit = 0.0; }
503  else if(val > CLHEP::pi) { polarAngleLimit = CLHEP::pi; }
504  else                     { polarAngleLimit = val; }
505}
506
507//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
508
509inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
510{
511  return stepLimit;
512}
513
514//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
515
516inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val) 
517{
518  stepLimit = val;
519  if(val == fMinimal) { facrange = 0.2; }
520}
521
522//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523
524#endif
Note: See TracBrowser for help on using the repository browser.