]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TGeant4/TG4G3CutVector.cxx
update to geant4 4.0 (not yet released)
[u/mrichter/AliRoot.git] / TGeant4 / TG4G3CutVector.cxx
CommitLineData
2817d3e2 1// $Id$
0375ab40 2// Category: global
2817d3e2 3//
e5967ab3 4// Author: I. Hrivnacova
5//
6// Class TG4G3CutVector
7// --------------------
2817d3e2 8// See the class description in the header file.
9
0375ab40 10#include "TG4G3CutVector.h"
2817d3e2 11#include "TG4G3Defaults.h"
12
13#include <G4Track.hh>
14#include <G4ParticleDefinition.hh>
15#include <G4VProcess.hh>
16
e5967ab3 17const G4double TG4G3CutVector::fgkDCUTEOff = 10. * TeV;
18const G4double TG4G3CutVector::fgkDCUTMOff = 10. * TeV;
19const G4double TG4G3CutVector::fgkTolerance = 1. * keV;
20TG4StringVector TG4G3CutVector::fgCutNameVector;
21
72095f7c 22//_____________________________________________________________________________
0375ab40 23TG4G3CutVector::TG4G3CutVector()
e5967ab3 24 : fDeltaRaysOn(true)
2817d3e2 25{
e5967ab3 26 // fill name vector
27 if (fgCutNameVector.size() == 0) FillCutNameVector();
28
2817d3e2 29 // initialize fCutVector
e5967ab3 30 for (G4int i=0; i<kTOFMAX; i++) fCutVector.push_back(0.);
31 fCutVector.push_back(DBL_MAX);
2817d3e2 32}
33
72095f7c 34//_____________________________________________________________________________
0375ab40 35TG4G3CutVector::TG4G3CutVector(const TG4G3CutVector& right)
e5967ab3 36 : fCutVector(right.fCutVector.size())
2817d3e2 37{
58c0119e 38 // copy stuff
39 *this = right;
2817d3e2 40}
41
72095f7c 42//_____________________________________________________________________________
0375ab40 43TG4G3CutVector::~TG4G3CutVector() {
2817d3e2 44//
2817d3e2 45}
46
47// operators
48
72095f7c 49//_____________________________________________________________________________
0375ab40 50TG4G3CutVector& TG4G3CutVector::operator=(const TG4G3CutVector& right)
2817d3e2 51{
52 // check assignement to self
53 if (this == &right) return *this;
54
e5967ab3 55 // copy fCutVector
56 for (G4int i=0; i<kNoG3Cuts; i++) fCutVector[i] = right.fCutVector[i];
57
58 fDeltaRaysOn = right.fDeltaRaysOn;
2817d3e2 59
60 return *this;
61}
62
72095f7c 63//_____________________________________________________________________________
0375ab40 64G4double TG4G3CutVector::operator[](G4int index) const
2817d3e2 65{
66//
67 if (index < kNoG3Cuts)
e5967ab3 68 return fCutVector[index];
2817d3e2 69 else {
70 TG4Globals::Exception(
0375ab40 71 "TG4G3CutVector::operator[]: index out of the vector scope");
2817d3e2 72 return 0.;
73 }
74}
75
e5967ab3 76// private methods
77
78//_____________________________________________________________________________
79void TG4G3CutVector::FillCutNameVector()
80{
81// Defines fCutNameVector.
82// ---
83
84 fgCutNameVector.push_back("CUTGAM");
85 fgCutNameVector.push_back("CUTELE");
86 fgCutNameVector.push_back("CUTNEU");
87 fgCutNameVector.push_back("CUTHAD");
88 fgCutNameVector.push_back("CUTMUO");
89 fgCutNameVector.push_back("BCUTE");
90 fgCutNameVector.push_back("BCUTM");
91 fgCutNameVector.push_back("DCUTE");
92 fgCutNameVector.push_back("DCUTM");
93 fgCutNameVector.push_back("PPCUTM");
94 fgCutNameVector.push_back("TOFMAX");
95 fgCutNameVector.push_back("NONE");
96}
97
2817d3e2 98// public methods
99
72095f7c 100//_____________________________________________________________________________
e5967ab3 101TG4G3Cut TG4G3CutVector::GetCut(const G4String& cutName)
102{
103// Retrieves corresponding TG4G3Cut constant from the cutName.
104// ---
105
106 if (cutName == fgCutNameVector[kCUTGAM]) return kCUTGAM;
107 else if (cutName == fgCutNameVector[kCUTELE]) return kCUTELE;
108 else if (cutName == fgCutNameVector[kCUTNEU]) return kCUTNEU;
109 else if (cutName == fgCutNameVector[kCUTHAD]) return kCUTHAD;
110 else if (cutName == fgCutNameVector[kCUTMUO]) return kCUTMUO;
111 else if (cutName == fgCutNameVector[kBCUTE]) return kBCUTE;
112 else if (cutName == fgCutNameVector[kBCUTM]) return kBCUTM;
113 else if (cutName == fgCutNameVector[kDCUTE]) return kDCUTE;
114 else if (cutName == fgCutNameVector[kDCUTM]) return kDCUTM;
115 else if (cutName == fgCutNameVector[kPPCUTM]) return kPPCUTM;
116 else if (cutName == fgCutNameVector[kTOFMAX]) return kTOFMAX;
117 else return kNoG3Cuts;
118}
119
120//_____________________________________________________________________________
121const G4String& TG4G3CutVector::GetCutName(TG4G3Cut cut)
122{
123// Returns name of a specified cut.
124// ---
125
126 // fill name vector
127 if (fgCutNameVector.size() == 0) TG4G3CutVector::FillCutNameVector();
128
129 return fgCutNameVector[cut];
130}
131
132//_____________________________________________________________________________
133void TG4G3CutVector::SetCut(TG4G3Cut cut, G4double cutValue)
2817d3e2 134{
135// Sets the cutValue for the specified cut.
136// ---
137
e5967ab3 138 if (cut>=kNoG3Cuts) {
2817d3e2 139 TG4Globals::Exception(
0375ab40 140 "TG4G3CutVector::SetG3Cut: Inconsistent cut.");
2817d3e2 141 }
e5967ab3 142
143 fCutVector[cut] = cutValue;
2817d3e2 144}
145
72095f7c 146//_____________________________________________________________________________
0375ab40 147void TG4G3CutVector::SetG3Defaults()
2817d3e2 148{
149// Sets G3 default values for all cuts.
150// ---
151
e5967ab3 152 for (G4int i=0; i<kNoG3Cuts; i++)
153 fCutVector[i] = TG4G3Defaults::Instance()->CutValue(i);
2817d3e2 154}
155
72095f7c 156//_____________________________________________________________________________
e5967ab3 157G4bool TG4G3CutVector::Update(const TG4G3CutVector& vector)
2817d3e2 158{
e5967ab3 159// Update only those values that have not yet been set ( =0.)
160// with given vector.
161// Returns true if some value was modified.
2817d3e2 162// ---
163
e5967ab3 164 G4bool result = false;
2817d3e2 165
e5967ab3 166 for (G4int i=0; i<kNoG3Cuts; i++)
167 if (fCutVector[i] < fgkTolerance ) {
168 fCutVector[i] = vector[i];
169 result = true;
170 }
171
172 return result;
2817d3e2 173}
174
e5967ab3 175//_____________________________________________________________________________
176void TG4G3CutVector::Print() const
177{
178// Prints the cuts.
179// ---
180
181 G4cout << " G3 cut vector:" << G4endl;
182 if (fDeltaRaysOn)
183 G4cout << " Delta rays On" << G4endl;
184 else
185 G4cout << " Delta rays Off" << G4endl;
186
187 // energy cuts
188 for (G4int i=0; i<kTOFMAX; i++) {
189
190 G4cout << " " << fgCutNameVector[i] << " cut value (MeV): ";
191
192 if (i == kDCUTE && !fDeltaRaysOn)
193 G4cout << fgkDCUTEOff/MeV << G4endl;
194 else if (i == kDCUTM && !fDeltaRaysOn)
195 G4cout << fgkDCUTMOff/MeV << G4endl;
196 else
197 G4cout << fCutVector[i]/MeV << G4endl;
198 }
199
200 // time cut
201 G4cout << " " << fgCutNameVector[kTOFMAX] << " cut value (s): "
202 << fCutVector[kTOFMAX]/s << G4endl;
203}
204
72095f7c 205//_____________________________________________________________________________
0375ab40 206G4double TG4G3CutVector::GetMinEkineForGamma(const G4Track& track) const
2817d3e2 207{
208// Returns the cut value for gamma.
209// (Cut is not applied for "opticalphoton"
210// as it is treated in G4 as a particle different
211// from "gamma" in G4.)
212// ---
213
2817d3e2 214 const G4VProcess* kpCreatorProcess = track.GetCreatorProcess();
215 G4String processName = "";
216 if (kpCreatorProcess) processName = kpCreatorProcess->GetProcessName();
217
218 if ((processName == "eBrem") || (processName == "IeBrem")) {
e5967ab3 219 return fCutVector[kBCUTE];
2817d3e2 220 }
e5967ab3 221 else if ((processName == "MuBrems") || (processName == "IMuBrems")) {
222 //(processName == "//hBrems")|| (processName == "//IhBrems")
2817d3e2 223 // hadron Brehmstrahlung is not defined in G4
e5967ab3 224 return fCutVector[kBCUTM];
2817d3e2 225 }
226 else {
e5967ab3 227 return fCutVector[kCUTGAM];
2817d3e2 228 }
229}
230
72095f7c 231//_____________________________________________________________________________
0375ab40 232G4double TG4G3CutVector::GetMinEkineForElectron(const G4Track& track) const
2817d3e2 233{
e5967ab3 234// Returns the cut value for e-.
2817d3e2 235// ---
236
237 const G4VProcess* kpCreatorProcess = track.GetCreatorProcess();
238 G4String processName = "";
239 if (kpCreatorProcess) processName = kpCreatorProcess->GetProcessName();
240
241 if ((processName == "eIoni") || (processName == "IeIoni")) {
e5967ab3 242 // delta rays by e-, e+
243
244 if (fDeltaRaysOn)
245 return fCutVector[kDCUTE];
246 else
247 return fgkDCUTEOff;
2817d3e2 248 }
e5967ab3 249 else if ((processName == "MuIoni") || (processName == "IMuIoni") ||
250 (processName == "hIoni") || (processName == "IhIoni")) {
251 // delta rays by other particles (mu, hadron)
252
253 if (fDeltaRaysOn)
254 return fCutVector[kDCUTM];
255 else
256 return fgkDCUTMOff;
2817d3e2 257 }
e5967ab3 258 else if (processName == "MuPairProd") {
259 //direct pair production by muons
260
261 return fCutVector[kPPCUTM];
262 }
2817d3e2 263 else {
e5967ab3 264 return fCutVector[kCUTELE];
2817d3e2 265 }
266}
267
72095f7c 268//_____________________________________________________________________________
e5967ab3 269G4double TG4G3CutVector::GetMinEkineForEplus(const G4Track& track) const
270{
271// Returns the cut value for e+.
272// ---
273
274 const G4VProcess* kpCreatorProcess = track.GetCreatorProcess();
275 G4String processName = "";
276 if (kpCreatorProcess) processName = kpCreatorProcess->GetProcessName();
277
278 if (processName == "MuPairProd") {
279 //direct pair production by muons
280 return fCutVector[kPPCUTM];
281 }
282 else
283 return fCutVector[kCUTELE];
284}
285
286//_____________________________________________________________________________
287G4double TG4G3CutVector::GetMinEkineForChargedHadron(const G4Track& track) const
2817d3e2 288{
289// Returns the cut value for charged hadron.
290// ---
291
e5967ab3 292 return fCutVector[kCUTHAD];
2817d3e2 293}
294
72095f7c 295//_____________________________________________________________________________
0375ab40 296G4double TG4G3CutVector::GetMinEkineForNeutralHadron(const G4Track& track) const
2817d3e2 297{
298// Returns the cut value for neutral hadron.
299// ---
300
e5967ab3 301 return fCutVector[kCUTNEU];
2817d3e2 302}
303
72095f7c 304//_____________________________________________________________________________
0375ab40 305G4double TG4G3CutVector::GetMinEkineForMuon(const G4Track& track) const
2817d3e2 306{
307// Returns the cut value for neutral muon.
308// ---
309
e5967ab3 310 return fCutVector[kCUTMUO];
2817d3e2 311}
312
72095f7c 313//_____________________________________________________________________________
0375ab40 314G4double TG4G3CutVector::GetMinEkineForOther(const G4Track& track) const
2817d3e2 315{
316// Returns 0.
317// ---
318
319 return 0.;
320}
e5967ab3 321//_____________________________________________________________________________
322G4bool TG4G3CutVector::IsCut() const
323{
324// Returns true if any of cuts is set.
325// ---
2817d3e2 326
e5967ab3 327 for (G4int i=0; i<kNoG3Cuts; i++)
328 if (fCutVector[i] > fgkTolerance ) return true;
329
330 return false;
331}