Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4Limits.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2014 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
15#include "TG4Limits.h"
16
17#include <G4SystemOfUnits.hh>
18
19G4int TG4Limits::fgCounter = 0;
20
21//_____________________________________________________________________________
23 const TG4G3CutVector& cuts, const TG4G3ControlVector& controls)
24 : G4UserLimits(),
25 // default values of G4UserLimits data members are set:
26 // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
27 // fMinEkine(0.), fMinRange(0.)
28 fName(""),
29 fIsCut(false),
30 fIsControl(false),
31 fCutVector(cuts),
32 fControlVector(),
33 fDefaultMaxStep(DBL_MAX)
34{
36
37 Initialize(cuts, controls);
38}
39
40//_____________________________________________________________________________
41TG4Limits::TG4Limits(const G4String& name, const TG4G3CutVector& cuts,
42 const TG4G3ControlVector& controls)
43 : G4UserLimits(),
44 // default values of G4UserLimits data members are set:
45 // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
46 // fMinEkine(0.), fMinRange(0.)
47 fName(name),
48 fIsCut(false),
49 fIsControl(false),
50 fCutVector(cuts),
51 fControlVector(),
52 fDefaultMaxStep(DBL_MAX)
53{
55
56 Initialize(cuts, controls);
57}
58
59//_____________________________________________________________________________
61 const TG4G3ControlVector& controls)
62 : G4UserLimits(g4Limits),
63 fName(""),
64 fIsCut(false),
65 fIsControl(false),
66 fCutVector(cuts),
67 fControlVector(),
68 fDefaultMaxStep(DBL_MAX)
69{
71
72 Initialize(cuts, controls);
73}
74
75//_____________________________________________________________________________
77 : G4UserLimits(),
78 // default values of G4UserLimits data members are set:
79 // fMaxStep (DBL_MAX), fMaxTrack(DBL_MAX),fMaxTime(DBL_MAX),
80 // fMinEkine(0.), fMinRange(0.)
81 fName(""),
82 fIsCut(false),
83 fIsControl(false),
84 fCutVector(),
85 fControlVector(),
86 fDefaultMaxStep(DBL_MAX)
87{
89
90 ++fgCounter;
91}
92
93//_____________________________________________________________________________
95 : G4UserLimits(right),
96 fName(right.fName),
97 fIsCut(right.fIsCut),
98 fIsControl(right.fIsControl),
99 fCutVector(right.fCutVector),
100 fControlVector(right.fControlVector),
101 fDefaultMaxStep(right.fDefaultMaxStep)
102{
104
105 ++fgCounter;
106}
107
108//_____________________________________________________________________________
113
114//
115// operators
116//
117
118//_____________________________________________________________________________
120{
122
123 // check assignement to self
124 if (this == &right) return *this;
125
126 // base class assignement
127 G4UserLimits::operator=(right);
128
129 fName = right.fName;
130 fIsCut = right.fIsCut;
131 fIsControl = right.fIsControl;
132 fCutVector = right.fCutVector;
134
135 return *this;
136}
137
138//
139// private methods
140//
141
142//_____________________________________________________________________________
144 const TG4G3CutVector& cuts, const TG4G3ControlVector& controls)
145{
147
148 fMaxTime = cuts[kTOFMAX];
149
150 fControlVector.Update(controls);
151 // only controls different from passed controls (default) are set
152
155
156 ++fgCounter;
157}
158
159//
160// public methods
161//
162
163//_____________________________________________________________________________
164G4double TG4Limits::GetUserMinEkine(const G4Track& /*track*/)
165{
169
170 return fMinEkine;
171}
172
173//_____________________________________________________________________________
174void TG4Limits::SetG3Cut(TG4G3Cut cut, G4double cutValue)
175{
177
178 fCutVector.SetCut(cut, cutValue);
179 fIsCut = true;
180
181 if (cut == kTOFMAX) fMaxTime = cutValue;
182}
183
184//_____________________________________________________________________________
186 TG4G3Control control, TG4G3ControlValue controlValue)
187{
189
190 G4bool result = fControlVector.SetControl(control, controlValue, fCutVector);
191
192 if (result) fIsControl = true;
193}
194
195//_____________________________________________________________________________
197{
199
201 fIsCut = true;
202}
203
204//_____________________________________________________________________________
206{
210
211 // G4bool result = fCutVector.Update(cutVector);
212 // if (result) fIsCut = true;
213
214 fControlVector.Update(controls);
216
217 return fIsControl;
218}
219
220//_____________________________________________________________________________
228
229//_____________________________________________________________________________
231{
233
234 // G4cout << "TG4Limits::SetCurrentMaxStep: in " << GetName()
235 // << " step " << step << G4endl;
236
237 fMaxStep = step;
238}
239
240//_____________________________________________________________________________
242{
245
246 // G4cout << "TG4Limits::SetDefaultMaxAllowedStep: in " << GetName()
247 // << " step " << fMaxStep << G4endl;
248
249 fDefaultMaxStep = fMaxStep;
250}
251
252//_____________________________________________________________________________
254{
257
258 // G4cout << "TG4Limits::SetMaxAllowedStepBack: in " << GetName()
259 // << " step " << fDefaultMaxStep << G4endl;
260
261 fMaxStep = fDefaultMaxStep;
262}
263
264//_____________________________________________________________________________
266{
268
269 G4cout << "\"" << fName << "\" limits:" << G4endl;
270 G4cout << " Max step length (mm): " << fMaxStep / mm << G4endl;
271 G4cout << " Max track length (mm): " << fMaxTrack / mm << G4endl;
272 G4cout << " Max time (s) " << fMaxTime / s << G4endl;
273 G4cout << " Min kin. energy (MeV) " << fMinEkine / MeV << G4endl;
274 G4cout << " Min range (mm): " << fMinRange / mm << G4endl;
275
276 if (!fIsCut) G4cout << " No special cuts. " << G4endl;
277 if (!fIsControl) G4cout << " No special controls. " << G4endl;
278
279 if (fIsCut) fCutVector.Print();
281}
282
283//_____________________________________________________________________________
284G4double TG4Limits::GetMinEkineForGamma(const G4Track& track) const
285{
287
288 if (fIsCut)
289 return fCutVector.GetMinEkineForGamma(track);
290 else
291 return fMinEkine;
292}
293
294//_____________________________________________________________________________
295G4double TG4Limits::GetMinEkineForElectron(const G4Track& track) const
296{
298
299 if (fIsCut)
301 else
302 return fMinEkine;
303}
304
305//_____________________________________________________________________________
306G4double TG4Limits::GetMinEkineForChargedHadron(const G4Track& track) const
307{
309
310 if (fIsCut)
312 else
313 return fMinEkine;
314}
315
316//_____________________________________________________________________________
317G4double TG4Limits::GetMinEkineForNeutralHadron(const G4Track& track) const
318{
320
321 if (fIsCut)
323 else
324 return fMinEkine;
325}
326
327//_____________________________________________________________________________
328G4double TG4Limits::GetMinEkineForMuon(const G4Track& track) const
329{
331
332 if (fIsCut)
333 return fCutVector.GetMinEkineForMuon(track);
334 else
335 return fMinEkine;
336}
337
338//_____________________________________________________________________________
340{
343
344 if (fIsControl)
345 return fControlVector.GetControlValue(process);
346 else
347 return kUnsetControlValue;
348}
Definition of the TG4Limits class.
Vector of control process values with convenient set/get methods.
G4bool Update(const TG4G3ControlVector &vector)
G4bool SetControl(TG4G3Control control, TG4G3ControlValue controlValue, TG4G3CutVector &cuts)
static TG4G3ControlValue GetControlValue(G4int value, TG4G3Control control)
Vector of kinetic energy cut values with convenient set/get methods.
G4double GetMinEkineForChargedHadron(const G4Track &track) const
G4double GetMinEkineForGamma(const G4Track &track) const
G4double GetMinEkineForMuon(const G4Track &track) const
G4double GetMinEkineForNeutralHadron(const G4Track &track) const
void Print() const
G4bool IsCut() const
G4double GetMinEkineForElectron(const G4Track &track) const
void SetCut(TG4G3Cut cut, G4double cutValue)
Extended G4UserLimits class.
Definition TG4Limits.h:38
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
void SetMaxAllowedStepBack()
TG4G3CutVector fCutVector
the vector of G3 cut values
Definition TG4Limits.h:98
TG4Limits & operator=(const TG4Limits &right)
virtual ~TG4Limits()
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
void Print() const
void SetCurrentMaxAllowedStep(G4double step)
G4bool fIsCut
true if any cut value is set
Definition TG4Limits.h:96
static G4int fgCounter
counter
Definition TG4Limits.h:92
G4double GetMinEkineForChargedHadron(const G4Track &track) const
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.
@ kUnsetControlValue
value not set
@ kTOFMAX
Definition TG4G3Cut.h:77