]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4Limits.h
added formatting in PrintStatistics()
[u/mrichter/AliRoot.git] / TGeant4 / TG4Limits.h
1 // $Id$
2 // Category: global
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class TG4Limits
7 // ---------------
8 // G4UserLimits derived class extended with the
9 // vectors of kinetic energy cuts and control process flags
10 // data members.
11
12 #ifndef TG4_LIMITS_H
13 #define TG4_LIMITS_H
14
15 #include "TG4Globals.h"
16 #include "TG4G3Cut.h"
17 #include "TG4G3Control.h"
18 #include "TG4G3CutVector.h"
19 #include "TG4G3ControlVector.h"
20
21 #include <G4UserLimits.hh>
22
23 class G4VProcess;
24
25 class TG4Limits: public G4UserLimits
26 {
27   public:
28     TG4Limits(const TG4G3CutVector& cuts, const TG4G3ControlVector& controls);
29     TG4Limits(const G4String& name,
30               const TG4G3CutVector& cuts, const TG4G3ControlVector& controls);
31     TG4Limits(const TG4Limits& right);
32     virtual ~TG4Limits();
33     
34     // operators
35     TG4Limits& operator=(const TG4Limits& right);
36     
37     // static methods
38     static G4int GetNofLimits();
39
40     // set methods
41     void SetName(const G4String& name);
42     void SetG3Cut(TG4G3Cut cut, G4double cutValue);
43     void SetG3Control(TG4G3Control control, TG4G3ControlValue controlValue);
44     G4bool Update(const TG4G3ControlVector& controls);
45     void SetG3DefaultCuts();
46     void SetG3DefaultControls();
47     
48     // methods
49     void Print() const;
50
51     // get methods
52     G4String GetName() const;
53     const TG4G3CutVector* GetCutVector() const;
54     const TG4G3ControlVector* GetControlVector() const;
55     G4bool  IsCut() const;
56     G4bool  IsControl() const;
57     virtual G4double GetUserMinEkine(const G4Track& track);
58     G4double GetMinEkineForGamma(const G4Track& track) const;
59     G4double GetMinEkineForElectron(const G4Track& track) const;
60     G4double GetMinEkineForEplus(const G4Track& track) const;
61     G4double GetMinEkineForChargedHadron(const G4Track& track) const;
62     G4double GetMinEkineForNeutralHadron(const G4Track& track) const;
63     G4double GetMinEkineForMuon(const G4Track& track) const;
64     G4double GetMinEkineForOther(const G4Track& track) const;
65     TG4G3ControlValue GetControl(G4VProcess* process) const; 
66
67   protected:
68     TG4Limits();
69
70   private:
71     // methods
72     void Initialize(const TG4G3CutVector& cuts, 
73                     const TG4G3ControlVector& controls);
74   
75     // static data members
76     static const G4double  fgkDefaultMaxStep; // default max step value
77     static G4int           fgCounter;         // counter 
78
79     // data members
80     G4String            fName;         //name
81     G4bool              fIsCut;        //true if any cut value is set
82     G4bool              fIsControl;    //true if any control value is set
83     TG4G3CutVector      fCutVector;    //TG4CutVector
84     TG4G3ControlVector  fControlVector;//TG4ControlVector
85 };
86
87 // inline methods
88
89 inline G4int TG4Limits::GetNofLimits()
90 { return fgCounter; }
91
92 inline G4bool TG4Limits::IsCut() const  
93 { return fIsCut; }
94
95 inline G4bool TG4Limits::IsControl() const 
96 { return fIsControl; }
97
98 inline void TG4Limits::SetName(const G4String& name) 
99 { fName = name; }
100
101 inline G4String TG4Limits::GetName() const
102 { return fName; }
103
104 inline const TG4G3CutVector* TG4Limits::GetCutVector() const
105 { return &fCutVector; }
106
107 inline const TG4G3ControlVector* TG4Limits::GetControlVector() const
108 { return &fControlVector; }
109
110 #endif //TG4_USER_LIMITS_H
111
112
113