Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4Limits.h
Go to the documentation of this file.
1#ifndef TG4_LIMITS_H
2#define TG4_LIMITS_H
3
4//------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2014 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18#include "TG4G3Control.h"
19#include "TG4G3ControlVector.h"
20#include "TG4G3Cut.h"
21#include "TG4G3CutVector.h"
22#include "TG4Globals.h"
23
24#include <G4UserLimits.hh>
25
26class G4VProcess;
27
36
37class TG4Limits : public G4UserLimits
38{
39 public:
40 TG4Limits(const TG4G3CutVector& cuts, const TG4G3ControlVector& controls);
41 TG4Limits(const G4String& name, const TG4G3CutVector& cuts,
42 const TG4G3ControlVector& controls);
43 TG4Limits(const G4UserLimits& g4Limits, const TG4G3CutVector& cuts,
44 const TG4G3ControlVector& controls);
45 TG4Limits(const TG4Limits& right);
46 virtual ~TG4Limits();
47
48 // operators
49 TG4Limits& operator=(const TG4Limits& right);
50
51 // static methods
52 static G4int GetNofLimits();
53
54 // set methods
55 void SetName(const G4String& name);
56 void SetG3Cut(TG4G3Cut cut, G4double cutValue);
57 void SetG3Control(TG4G3Control control, TG4G3ControlValue controlValue);
58 G4bool Update(const TG4G3ControlVector& controls);
59 void SetG3DefaultCuts();
61 void SetCurrentMaxAllowedStep(G4double step);
64
65 // methods
66 void Print() const;
67
68 // get methods
69 G4String GetName() const;
70 G4double GetMaxUserStep() const;
71 const TG4G3CutVector* GetCutVector() const;
73 G4bool IsCut() const;
74 G4bool IsControl() const;
75 virtual G4double GetUserMinEkine(const G4Track& track);
76 G4double GetMinEkineForGamma(const G4Track& track) const;
77 G4double GetMinEkineForElectron(const G4Track& track) const;
78 G4double GetMinEkineForChargedHadron(const G4Track& track) const;
79 G4double GetMinEkineForNeutralHadron(const G4Track& track) const;
80 G4double GetMinEkineForMuon(const G4Track& track) const;
82
83 private:
85 TG4Limits();
86
87 // methods
88 void Initialize(
89 const TG4G3CutVector& cuts, const TG4G3ControlVector& controls);
90
91 // static data members
92 static G4int fgCounter;
93
94 // data members
95 G4String fName;
96 G4bool fIsCut;
97 G4bool fIsControl;
101};
102
103// inline methods
104
106{
108 return fgCounter;
109}
110
111inline G4bool TG4Limits::IsCut() const
112{
114 return fIsCut;
115}
116
117inline G4bool TG4Limits::IsControl() const
118{
120 return fIsControl;
121}
122
123inline void TG4Limits::SetName(const G4String& name)
124{
126 fName = name;
127}
128
129inline G4String TG4Limits::GetName() const
130{
132 return fName;
133}
134
135inline G4double TG4Limits::GetMaxUserStep() const
136{
138 return fMaxStep;
139}
140
142{
144 return &fCutVector;
145}
146
148{
150 return &fControlVector;
151}
152
153#endif // TG4_USER_LIMITS_H
Definition of the TG4G3ControlVector class.
Definition of the enumerations TG4G3Control, TG4G3ControlValue.
Definition of the TG4G3CutVector class.
Definition of the enumeration TG4G3Cut.
Definition of the TG4Globals class and basic container types.
Vector of control process values with convenient set/get methods.
Vector of kinetic energy cut values with convenient set/get methods.
Extended G4UserLimits class.
Definition TG4Limits.h:38
G4String GetName() const
Definition TG4Limits.h:129
G4bool Update(const TG4G3ControlVector &controls)
void SetDefaultMaxAllowedStep()
void SetG3Cut(TG4G3Cut cut, G4double cutValue)
TG4G3ControlValue GetControl(G4VProcess *process) const
void SetG3Control(TG4G3Control control, TG4G3ControlValue controlValue)
virtual G4double GetUserMinEkine(const G4Track &track)
TG4G3ControlVector fControlVector
the vector of G3 control values
Definition TG4Limits.h:99
G4double GetMinEkineForNeutralHadron(const G4Track &track) const
G4double GetMinEkineForMuon(const G4Track &track) const
G4double GetMaxUserStep() const
Definition TG4Limits.h:135
void SetMaxAllowedStepBack()
TG4G3CutVector fCutVector
the vector of G3 cut values
Definition TG4Limits.h:98
TG4Limits & operator=(const TG4Limits &right)
virtual ~TG4Limits()
const TG4G3CutVector * GetCutVector() const
Definition TG4Limits.h:141
G4double GetMinEkineForElectron(const G4Track &track) const
TG4Limits()
Not implemented.
Definition TG4Limits.cxx:76
void SetG3DefaultCuts()
void SetG3DefaultControls()
G4String fName
name
Definition TG4Limits.h:95
G4bool fIsControl
true if any control value is set
Definition TG4Limits.h:97
void Initialize(const TG4G3CutVector &cuts, const TG4G3ControlVector &controls)
G4double fDefaultMaxStep
the default max step value
Definition TG4Limits.h:100
G4double GetMinEkineForGamma(const G4Track &track) const
G4bool IsCut() const
Definition TG4Limits.h:111
static G4int GetNofLimits()
Definition TG4Limits.h:105
void SetName(const G4String &name)
Definition TG4Limits.h:123
void Print() const
void SetCurrentMaxAllowedStep(G4double step)
G4bool fIsCut
true if any cut value is set
Definition TG4Limits.h:96
G4bool IsControl() const
Definition TG4Limits.h:117
static G4int fgCounter
counter
Definition TG4Limits.h:92
G4double GetMinEkineForChargedHadron(const G4Track &track) const
const TG4G3ControlVector * GetControlVector() const
Definition TG4Limits.h:147
TG4G3ControlValue
Enumeration for G3 processes control values.
TG4G3Cut
Enumeration for G3 types of kinetic energy cuts.
Definition TG4G3Cut.h:22
TG4G3Control
Enumeration for G3 types of physics processes controls.