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

Last change on this file since 1055 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

File size: 17.9 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.56 2009/04/07 18:39:47 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-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  G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
130
131  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
132
133  // Store PhysicsTable in a file.
134  // Return false in case of failure at I/O
135  G4bool StorePhysicsTable(const G4ParticleDefinition*,
136                           const G4String& directory,
137                           G4bool ascii = false);
138
139  // Retrieve Physics from a file.
140  // (return true if the Physics Table can be build by using file)
141  // (return false if the process has no functionality or in case of failure)
142  // File name should is constructed as processName+particleName and the
143  // should be placed under the directory specifed by the argument.
144  G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
145                              const G4String& directory,
146                              G4bool ascii);
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.
151  G4double AlongStepGetPhysicalInteractionLength(
152                                            const G4Track&,
153                                            G4double  previousStepSize,
154                                            G4double  currentMinimalStep,
155                                            G4double& currentSafety,
156                                            G4GPILSelection* selection);
157
158  // The function overloads the corresponding function of the base
159  // class.
160  G4double PostStepGetPhysicalInteractionLength(
161                                            const G4Track&,
162                                            G4double  previousStepSize,
163                                            G4ForceCondition* condition);
164
165  // This method does not used for tracking, it is intended only for tests
166  inline G4double ContinuousStepLimit(const G4Track& track,
167                                      G4double previousStepSize,
168                                      G4double currentMinimalStep,
169                                      G4double& currentSafety);
170
171  //------------------------------------------------------------------------
172  // Specific methods to build and access Physics Tables
173  //------------------------------------------------------------------------
174
175  // Build empty Physics Vector
176  G4PhysicsVector* PhysicsVector(const G4MaterialCutsCouple*);
177
178  inline void SetBinning(G4int nbins);
179  inline G4int Binning() const;
180
181  inline void SetMinKinEnergy(G4double e);
182  inline G4double MinKinEnergy() const;
183
184  inline void SetMaxKinEnergy(G4double e);
185  inline G4double MaxKinEnergy() const;
186
187  inline void SetBuildLambdaTable(G4bool val);
188
189  inline G4PhysicsTable* LambdaTable() const;
190
191  // access particle type
192  inline const G4ParticleDefinition* Particle() const;
193
194  //------------------------------------------------------------------------
195  // Specific methods to set, access, modify models
196  //------------------------------------------------------------------------
197
198protected:
199  // Select model in run time
200  inline G4VEmModel* SelectModel(G4double kinEnergy);
201
202public:
203  // Select model in run time
204  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 
205                                            size_t& idxRegion) const;
206
207  // Add model for region, smaller value of order defines which
208  // model will be selected for a given energy interval 
209  void AddEmModel(G4int order, G4VEmModel*, const G4Region* region = 0);
210
211  // Access to models by index
212  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
213
214  //------------------------------------------------------------------------
215  // Get/Set parameters for simulation of multiple scattering
216  //------------------------------------------------------------------------
217
218  inline G4bool LateralDisplasmentFlag() const;
219  inline void SetLateralDisplasmentFlag(G4bool val);
220
221  inline G4double Skin() const;
222  inline void SetSkin(G4double val);
223
224  inline G4double RangeFactor() const;
225  inline void SetRangeFactor(G4double val);
226
227  inline G4double GeomFactor() const;
228  inline void SetGeomFactor(G4double val);
229
230  inline G4double PolarAngleLimit() const;
231  inline void SetPolarAngleLimit(G4double val);
232
233  inline G4MscStepLimitType StepLimitType() const;
234  inline void SetStepLimitType(G4MscStepLimitType val);
235
236  //------------------------------------------------------------------------
237  // Run time methods
238  //------------------------------------------------------------------------
239
240protected:
241
242  // This method is not used for tracking, it returns mean free path value
243  G4double GetMeanFreePath(const G4Track& track,
244                           G4double,
245                           G4ForceCondition* condition);
246
247  // This method is not used for tracking, it returns step limit
248  G4double GetContinuousStepLimit(const G4Track& track,
249                                  G4double previousStepSize,
250                                  G4double currentMinimalStep,
251                                  G4double& currentSafety);
252
253  // This method returns inversed transport cross section
254  inline G4double GetLambda(const G4ParticleDefinition* p, 
255                            G4double& kineticEnergy);
256
257  // This method is used for tracking, it returns step limit
258  inline G4double GetMscContinuousStepLimit(const G4Track& track,
259                                            G4double scaledKinEnergy,
260                                            G4double currentMinimalStep,
261                                            G4double& currentSafety);
262
263  // defines current material in run time
264  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
265
266  inline const G4MaterialCutsCouple* CurrentMaterialCutsCouple() const; 
267
268private:
269
270  // hide  assignment operator
271  G4VMultipleScattering(G4VMultipleScattering &);
272  G4VMultipleScattering & operator=(const G4VMultipleScattering &right);
273
274  // ======== Parameters of the class fixed at construction =========
275
276  G4EmModelManager*           modelManager;
277  G4bool                      buildLambdaTable;
278
279  // ======== Parameters of the class fixed at initialisation =======
280
281  G4PhysicsTable*             theLambdaTable;
282  const G4ParticleDefinition* firstParticle;
283
284  G4MscStepLimitType          stepLimit;
285
286  G4double                    minKinEnergy;
287  G4double                    maxKinEnergy;
288  G4double                    skin;
289  G4double                    facrange;
290  G4double                    facgeom;
291  G4double                    polarAngleLimit;
292
293  G4int                       nBins;
294
295  G4bool                      latDisplasment;
296
297  // ======== Cashed values - may be state dependent ================
298
299protected:
300
301  G4GPILSelection             valueGPILSelectionMSC;
302  G4ParticleChangeForMSC      fParticleChange;
303
304private:
305
306  G4VMscModel*                currentModel;
307
308  // cache
309  const G4ParticleDefinition* currentParticle;
310  const G4MaterialCutsCouple* currentCouple;
311  size_t                      currentMaterialIndex;
312
313};
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
317
318inline G4double G4VMultipleScattering::ContinuousStepLimit(
319                                       const G4Track& track,
320                                       G4double previousStepSize,
321                                       G4double currentMinimalStep,
322                                       G4double& currentSafety)
323{
324  return GetMscContinuousStepLimit(track,previousStepSize,currentMinimalStep,
325                                   currentSafety);
326}
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329
330inline void G4VMultipleScattering::SetBinning(G4int nbins)
331{
332  nBins = nbins;
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
336
337inline G4int G4VMultipleScattering::Binning() const
338{
339  return nBins;
340}
341
342//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343
344inline void G4VMultipleScattering::SetMinKinEnergy(G4double e)
345{
346  minKinEnergy = e;
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350
351inline G4double G4VMultipleScattering::MinKinEnergy() const
352{
353  return minKinEnergy;
354}
355
356//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
357
358inline void G4VMultipleScattering::SetMaxKinEnergy(G4double e)
359{
360  maxKinEnergy = e;
361}
362
363//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
364
365inline G4double G4VMultipleScattering::MaxKinEnergy() const
366{
367  return maxKinEnergy;
368}
369
370//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
371
372inline  void G4VMultipleScattering::SetBuildLambdaTable(G4bool val)
373{
374  buildLambdaTable = val;
375}
376
377//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
378
379inline G4PhysicsTable* G4VMultipleScattering::LambdaTable() const
380{
381  return theLambdaTable;
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
385
386inline  const G4ParticleDefinition* G4VMultipleScattering::Particle() const
387{
388  return currentParticle;
389}
390
391//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
392
393inline G4VEmModel* G4VMultipleScattering::SelectModel(G4double kinEnergy)
394{
395  return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
396}
397
398//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
399
400inline G4VEmModel* G4VMultipleScattering::SelectModelForMaterial(
401                   G4double kinEnergy, size_t& idxRegion) const
402{
403  return modelManager->SelectModel(kinEnergy, idxRegion);
404}
405
406//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
407
408inline  G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
409{
410  return latDisplasment;
411}
412
413//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
414
415inline  void G4VMultipleScattering::SetLateralDisplasmentFlag(G4bool val)
416{
417  latDisplasment = val;
418}
419
420//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
421
422inline  G4double G4VMultipleScattering::Skin() const
423{
424  return skin;
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
429inline  void G4VMultipleScattering::SetSkin(G4double val)
430{
431  if(val < 1.0) skin = 0.0;
432  else          skin = val;
433}
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
436
437inline  G4double G4VMultipleScattering::RangeFactor() const
438{
439  return facrange;
440}
441
442//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
443
444inline  void G4VMultipleScattering::SetRangeFactor(G4double val)
445{
446  if(val > 0.0) facrange = val;
447}
448
449//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
450
451inline  G4double G4VMultipleScattering::GeomFactor() const
452{
453  return facgeom;
454}
455
456//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
457
458inline  void G4VMultipleScattering::SetGeomFactor(G4double val)
459{
460  if(val > 0.0) facgeom = val;
461}
462
463//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
464
465inline  G4double G4VMultipleScattering::PolarAngleLimit() const
466{
467  return polarAngleLimit;
468}
469
470//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
471
472inline  void G4VMultipleScattering::SetPolarAngleLimit(G4double val)
473{
474  if(val < 0.0)     polarAngleLimit = 0.0;
475  else if(val > pi) polarAngleLimit = pi;
476  else              polarAngleLimit = val;
477}
478
479//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480
481inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
482{
483  return stepLimit;
484}
485
486//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
487
488inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val) 
489{
490  stepLimit = val;
491  if(val == fMinimal) facrange = 0.2;
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
495
496inline 
497G4double G4VMultipleScattering::GetLambda(const G4ParticleDefinition* p, 
498                                          G4double& e)
499{
500  G4double x;
501  if(theLambdaTable) {
502    G4bool b;
503    x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
504  } else {
505    x = currentModel->CrossSection(currentCouple,p,e);
506  }
507  if(x > DBL_MIN) x = 1./x;
508  else            x = DBL_MAX; 
509  return x;
510}
511
512//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
513
514inline G4double G4VMultipleScattering::GetMscContinuousStepLimit(
515                                          const G4Track& track,
516                                          G4double scaledKinEnergy,
517                                          G4double currentMinimalStep,
518                                          G4double&)
519{
520  G4double x = currentMinimalStep;
521  DefineMaterial(track.GetMaterialCutsCouple());
522  currentModel = static_cast<G4VMscModel*>(SelectModel(scaledKinEnergy));
523  if(x > 0.0 && scaledKinEnergy > 0.0) {
524    G4double tPathLength = 
525      currentModel->ComputeTruePathLengthLimit(track, theLambdaTable, x);
526    if (tPathLength < x) valueGPILSelectionMSC = CandidateForSelection;
527    x = currentModel->ComputeGeomPathLength(tPathLength); 
528    //  G4cout << "tPathLength= " << tPathLength
529    //         << " stepLimit= " << x
530    //        << " currentMinimalStep= " << currentMinimalStep<< G4endl;
531  }
532  return x;
533}
534
535//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
536
537inline 
538void G4VMultipleScattering::DefineMaterial(const G4MaterialCutsCouple* couple)
539{
540  if(couple != currentCouple) {
541    currentCouple   = couple;
542    currentMaterialIndex = couple->GetIndex();
543  }
544}
545
546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547
548inline const G4MaterialCutsCouple* 
549G4VMultipleScattering::CurrentMaterialCutsCouple() const
550{
551  return currentCouple;
552} 
553
554//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
555
556#endif
Note: See TracBrowser for help on using the repository browser.