]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4Limits.cxx
Functions renamed to get a prefix PHOS
[u/mrichter/AliRoot.git] / TGeant4 / TG4Limits.cxx
1 // $Id$
2 // Category: global
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class TG4Limits
7 // ---------------
8 // See the class description in the header file.
9
10 #include "TG4Limits.h"
11
12 #include <globals.hh>
13
14 const G4double TG4Limits::fgkDefaultMaxStep = DBL_MAX;
15 G4int          TG4Limits::fgCounter = 0;
16
17 //_____________________________________________________________________________
18 TG4Limits::TG4Limits(const TG4G3CutVector& cuts, 
19                      const TG4G3ControlVector& controls)
20   : G4UserLimits(),              
21     // default values of G4UserLimits data members are set: 
22     // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
23     // fMinEkine(0.), fMinRange(0.)
24     fName(""),
25     fCutVector(cuts),
26     fControlVector(),
27     fIsCut(false),
28     fIsControl(false) 
29 {
30 //
31   Initialize(cuts, controls);
32 }
33
34 //_____________________________________________________________________________
35 TG4Limits::TG4Limits(const G4String& name, 
36                      const TG4G3CutVector& cuts, 
37                      const TG4G3ControlVector& controls)
38   : G4UserLimits(),              
39     // default values of G4UserLimits data members are set: 
40     // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
41     // fMinEkine(0.), fMinRange(0.)
42     fName(name),
43     fCutVector(cuts),
44     fControlVector(),
45     fIsCut(false),
46     fIsControl(false) 
47 {
48 //
49   Initialize(cuts, controls);
50 }
51
52 //_____________________________________________________________________________
53 TG4Limits::TG4Limits()
54   : G4UserLimits(),              
55     // default values of G4UserLimits data members are set: 
56     // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
57     // fMinEkine(0.), fMinRange(0.)
58     fName(""),
59     fIsCut(false),
60     fIsControl(false) 
61 {
62 //
63   fMaxStep = fgkDefaultMaxStep;
64   
65   ++fgCounter;
66 }
67
68 //_____________________________________________________________________________
69 TG4Limits::TG4Limits(const TG4Limits& right)
70   : G4UserLimits(right) 
71 {
72 //    
73   // copy stuff
74   *this = right;
75
76   ++fgCounter;
77 }  
78
79 //_____________________________________________________________________________
80 TG4Limits::~TG4Limits() {
81 //
82 }
83
84 // operators
85
86 //_____________________________________________________________________________
87 TG4Limits& TG4Limits::operator=(const TG4Limits& right)
88 {
89   // check assignement to self
90   if (this == &right) return *this;
91
92   // base class assignement
93   G4UserLimits::operator=(right);
94
95   fName = right.fName;
96   fIsCut = right.fIsCut;
97   fIsControl = right.fIsControl;
98   fCutVector  = right.fCutVector;
99   fControlVector = right.fControlVector;
100
101   return *this;  
102 }    
103           
104 // private methods
105
106 //_____________________________________________________________________________
107 void TG4Limits::Initialize(const TG4G3CutVector& cuts, 
108                            const TG4G3ControlVector& controls)
109 {                        
110 // Initialization.
111 // ---
112
113   fMaxStep = fgkDefaultMaxStep;
114   fMaxTime = cuts[kTOFMAX];
115
116   fControlVector.Update(controls);
117        // only controls different from passed controls (default) are set
118        
119   fIsCut = fCutVector.IsCut();
120   fIsControl = fControlVector.IsControl();     
121   
122   ++fgCounter;
123 }  
124   
125
126 // public methods
127
128 //_____________________________________________________________________________
129 G4double TG4Limits::GetUserMinEkine(const G4Track& track)
130 {
131 // Returns the kinetic energy cut parameter.
132 // !! The cuts values defined if fCutVector are applied
133 // only via TG4SpecialCuts process.
134 // ---
135
136   return fMinEkine;
137 }
138
139 //_____________________________________________________________________________
140 void TG4Limits::SetG3Cut(TG4G3Cut cut, G4double cutValue)
141 {
142 // Sets the cut value for the specified cut.
143 // ---
144
145   fCutVector.SetCut(cut, cutValue);
146   fIsCut = true;
147   
148   if (cut == kTOFMAX) fMaxTime = cutValue;
149 }
150
151 //_____________________________________________________________________________
152 void TG4Limits::SetG3Control(TG4G3Control control, 
153                              TG4G3ControlValue controlValue)
154 {
155 // Sets the process control value for the specified control.
156 // ---
157
158   G4bool result
159     = fControlVector.SetControl(control, controlValue, fCutVector);
160
161   if (result) fIsControl = true;
162   
163 }
164
165 //_____________________________________________________________________________
166 void TG4Limits::SetG3DefaultCuts()
167 {
168 // Sets the G3 default cut values for all cuts.
169 // ---
170
171   fCutVector.SetG3Defaults();
172   fIsCut = true;  
173 }
174
175 //_____________________________________________________________________________
176 G4bool TG4Limits::Update(const TG4G3ControlVector& controls)
177 {                      
178 // Updates controls in a special way.
179 // Returns true if some control in fControlVector is after update
180 // still set.
181 // ---
182
183   //G4bool result = fCutVector.Update(cutVector);
184   //if (result) fIsCut = true;
185
186   fControlVector.Update(controls);
187   fIsControl = fControlVector.IsControl();
188   
189   return fIsControl;
190 }
191
192 //_____________________________________________________________________________
193 void TG4Limits::SetG3DefaultControls()
194 {
195 // Sets the G3 default process control values for all flags.
196 // ---
197
198   fControlVector.SetG3Defaults();
199   fIsControl = true;
200 }
201
202 //_____________________________________________________________________________
203 void TG4Limits::Print() const
204 {
205 // Prints limits.
206 // ---
207
208    G4cout << "\"" << fName << "\"  limits:"<< G4endl;
209    G4cout << "  Max step length (mm):  " << fMaxStep/mm   << G4endl;
210    G4cout << "  Max track length (mm): " << fMaxTrack/mm  << G4endl;
211    G4cout << "  Max time (s)           " << fMaxTime/s    << G4endl;
212    G4cout << "  Min kin. energy (MeV)  " << fMinEkine/MeV << G4endl;
213    G4cout << "  Min range (mm):        " << fMinRange/mm  << G4endl;
214
215    if (!fIsCut)     G4cout << "  No special cuts. "<< G4endl;
216    if (!fIsControl) G4cout << "  No special controls. "<< G4endl;
217
218    if (fIsCut)     fCutVector.Print();
219    if (fIsControl) fControlVector.Print();
220 }
221
222 //_____________________________________________________________________________
223 G4double TG4Limits::GetMinEkineForGamma(const G4Track& track) const
224 {
225 // Returns the cut value for gamma.
226 // ---
227
228   if (fIsCut)
229     return fCutVector.GetMinEkineForGamma(track);
230   else 
231     return fMinEkine;
232 }
233
234 //_____________________________________________________________________________
235 G4double TG4Limits::GetMinEkineForElectron(const G4Track& track) const
236 {
237 // Returns the cut value for e-.
238 // ---
239
240   if (fIsCut)
241     return fCutVector.GetMinEkineForElectron(track);
242   else 
243     return fMinEkine;
244 }
245
246 //_____________________________________________________________________________
247 G4double TG4Limits::GetMinEkineForEplus(const G4Track& track) const
248 {
249 // Returns the cut value for e-.
250 // ---
251
252   if (fIsCut)
253     return fCutVector.GetMinEkineForEplus(track);
254   else 
255     return fMinEkine;
256 }
257
258 //_____________________________________________________________________________
259 G4double TG4Limits::GetMinEkineForChargedHadron(const G4Track& track) const
260 {
261 // Returns the cut value for charged hadron.
262 // ---
263
264   if (fIsCut)
265     return fCutVector.GetMinEkineForChargedHadron(track);
266   else 
267     return fMinEkine;
268 }
269
270 //_____________________________________________________________________________
271 G4double TG4Limits::GetMinEkineForNeutralHadron(const G4Track& track) const
272 {
273 // Returns the cut value for neutral hadron.
274 // ---
275
276   if (fIsCut)
277     return fCutVector.GetMinEkineForNeutralHadron(track);
278   else 
279     return fMinEkine;
280 }
281
282 //_____________________________________________________________________________
283 G4double TG4Limits::GetMinEkineForMuon(const G4Track& track) const
284 {
285 // Returns the cut value for neutral muon.
286 // ---
287
288   if (fIsCut)
289     return fCutVector.GetMinEkineForMuon(track);
290   else 
291     return fMinEkine;
292 }
293
294 //_____________________________________________________________________________
295 G4double TG4Limits::GetMinEkineForOther(const G4Track& track) const
296 {
297  // Returns 0.
298 // ---
299
300   if (fIsCut)
301     return fCutVector.GetMinEkineForOther(track);
302   else 
303     return fMinEkine;
304 }
305
306 //_____________________________________________________________________________
307 TG4G3ControlValue TG4Limits::GetControl(G4VProcess* process) const 
308 {
309 // Returns the flag value for the particle associated with
310 // the specified process.
311 // ---
312
313   if (fIsControl)
314     return fControlVector.GetControlValue(process);
315   else 
316     return kUnset;
317 }