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

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

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