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

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

update processes

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