Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4G3PhysicsManager.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 "TG4G3PhysicsManager.h"
16#include "TG4G3ControlVector.h"
17#include "TG4G3CutVector.h"
18#include "TG4G3Defaults.h"
19#include "TG4G3Units.h"
20
21#include <G4ParticleDefinition.hh>
22#include <G4ProcessTable.hh>
23#include <G4UImessenger.hh>
24#include <G4VProcess.hh>
25
27
28//_____________________________________________________________________________
30 : fCutVector(0),
31 fControlVector(0),
32 fIsCutVector(0),
33 fIsControlVector(0),
34 fG3Defaults(),
35 fLock(false)
36{
38
39 if (fgInstance) {
40 TG4Globals::Exception("TG4G3PhysicsManager", "TG4G3PhysicsManager",
41 "Cannot create two instances of singleton.");
42 }
43
44 fgInstance = this;
45
46 // create vector of kinetic energy cut values
48
49 // create vector of control process control values
51
52 // initialize fIsCutVector
54 G4int i;
55 for (i = 0; i < kNofParticlesWSP; i++) fIsCutVector->push_back(false);
56
57 // initialize fIsControlVector
59 for (i = 0; i < kNofParticlesWSP; i++) fIsControlVector->push_back(false);
60}
61
62//_____________________________________________________________________________
64{
66
67 delete fCutVector;
68 delete fControlVector;
69 delete fIsCutVector;
70 delete fIsControlVector;
71
72 fgInstance = 0;
73}
74
75//
76// private methods
77//
78
79//_____________________________________________________________________________
80void TG4G3PhysicsManager::SetCut(TG4G3Cut cut, G4double cutValue)
81{
83
84 // Do nothing if cutValue is not valid
85 if (!TG4G3CutVector::CheckCutValue(cut, cutValue)) return;
86
87 fCutVector->SetCut(cut, cutValue);
89}
90
91//_____________________________________________________________________________
93 TG4G3Control control, TG4G3ControlValue controlValue)
94{
96
97 if (control == kDRAY || control == kG3LOSS) {
100 }
101
102 fControlVector->SetControl(control, controlValue, *fCutVector);
103}
104
105//_____________________________________________________________________________
107{
109
110 switch (cut) {
111 case kCUTGAM:
112 (*fIsCutVector)[kGamma] = true;
113 break;
114 case kBCUTE:
115 (*fIsCutVector)[kGamma] = true;
116 break;
117 case kBCUTM:
118 (*fIsCutVector)[kGamma] = true;
119 break;
120 case kCUTELE:
121 (*fIsCutVector)[kElectron] = true;
122 break;
123 case kDCUTE:
124 (*fIsCutVector)[kElectron] = true;
125 break;
126 case kDCUTM:
127 (*fIsCutVector)[kElectron] = true;
128 break;
129 case kPPCUTM:
130 // (*fIsCutVector)[kElectron] = true;
131 // (*fIsCutVector)[kEplus] = true;
132 // PPCUTM is not applied with special cuts process
133 break;
134 case kCUTNEU:
135 (*fIsCutVector)[kNeutralHadron] = true;
136 break;
137 case kCUTHAD:
138 (*fIsCutVector)[kChargedHadron] = true;
139 break;
140 case kCUTMUO:
141 (*fIsCutVector)[kMuon] = true;
142 break;
143 default:
144 break;
145 }
146}
147
148//_____________________________________________________________________________
150{
153
154 switch (control) {
155 case kPAIR:
156 // gamma
157 (*fIsControlVector)[kGamma] = true;
158 break;
159 case kCOMP:
160 // gamma
161 (*fIsControlVector)[kGamma] = true;
162 break;
163 case kPHOT:
164 // gamma
165 (*fIsControlVector)[kGamma] = true;
166 break;
167 case kPFIS:
168 // gamma
169 (*fIsControlVector)[kGamma] = true;
170 break;
171 case kDRAY:
172 // all charged particles
173 (*fIsControlVector)[kElectron] = true;
174 (*fIsControlVector)[kEplus] = true;
175 (*fIsControlVector)[kChargedHadron] = true;
176 (*fIsControlVector)[kMuon] = true;
177 break;
178 case kANNI:
179 // e+ only
180 (*fIsControlVector)[kEplus] = true;
181 break;
182 case kBREM:
183 // e-/e+, muons
184 (*fIsControlVector)[kElectron] = true;
185 (*fIsControlVector)[kEplus] = true;
186 (*fIsControlVector)[kMuon] = true;
187 break;
188 case kHADR:
189 // hadrons
190 (*fIsControlVector)[kNeutralHadron] = true;
191 (*fIsControlVector)[kChargedHadron] = true;
192 break;
193 case kMUNU:
194 // muons
195 (*fIsControlVector)[kMuon] = true;
196 break;
197 case kDCAY:
198 // any
199 (*fIsControlVector)[kAny] = true;
200 break;
201 case kG3LOSS:
202 // all charged particles
203 (*fIsControlVector)[kElectron] = true;
204 (*fIsControlVector)[kEplus] = true;
205 (*fIsControlVector)[kChargedHadron] = true;
206 (*fIsControlVector)[kMuon] = true;
207 break;
208 case kMULS:
209 // all charged particles
210 (*fIsControlVector)[kElectron] = true;
211 (*fIsControlVector)[kEplus] = true;
212 (*fIsControlVector)[kChargedHadron] = true;
213 (*fIsControlVector)[kMuon] = true;
214 break;
215 default:
216 break;
217 }
218}
219
220//
221// public methods
222//
223
224//_____________________________________________________________________________
226{
229
230 if (fLock) {
231 TG4Globals::Exception("TG4PhysicsManager", "CheckLock",
232 "It is too late to change physics setup." + TG4Globals::Endl() +
233 "PhysicsManager has been already locked.");
234 }
235}
236
237//_____________________________________________________________________________
239 G4String name, G4double value, TG4G3Cut& cut)
240{
246
247 // convert cut name -> TG4G3Cut
248 cut = TG4G3CutVector::GetCut(name);
249
250 // add units
251 if (cut == kTOFMAX)
252 value *= TG4G3Units::Time();
253 else
254 value *= TG4G3Units::Energy();
255
256 // check cut value validity
257 if (!TG4G3CutVector::CheckCutValue(cut, value)) return false;
258
259 // set switch vector element only if the value
260 // is different from the value in cutVector and is valid
261 if (cut != kNoG3Cuts) {
262 // G4cout << "Comparing: "
263 // << value << " "
264 // << (*fCutVector)[cut] << " "
265 // << std::abs(value - (*fCutVector)[cut]) << " > "
266 // << TG4G3CutVector::Tolerance() << G4endl;
267
268 if (std::abs(value - (*fCutVector)[cut]) > TG4G3CutVector::Tolerance()) {
270 return true;
271 }
272 else
273 return false;
274 }
275 return false;
276}
277
278//_____________________________________________________________________________
280 G4double value, TG4G3Control& control, TG4G3ControlValue& controlValue)
281{
287
288 // convert control name -> TG4G3Control
289 control = TG4G3ControlVector::GetControl(name);
290
291 // convert double value -> TG4G3ControlValue
292 controlValue = TG4G3ControlVector::GetControlValue(value, control);
293
294 // set switch vector element only if the value
295 // is different from the value in controlVector
296 if (control != kNoG3Controls) {
297 if (controlValue != (*fControlVector)[control]) {
298 SwitchIsControlVector(control);
299 return true;
300 }
301 else
302 return false;
303 }
304 return false;
305}
306
307//_____________________________________________________________________________
309{
311
312 CheckLock();
314
315 for (G4int i = 0; i < kNofParticlesWSP; i++) (*fIsCutVector)[i] = true;
316}
317
318//_____________________________________________________________________________
320{
322
323 CheckLock();
325
326 for (G4int i = 0; i < kNofParticlesWSP; i++) (*fIsControlVector)[i] = true;
327}
328
329//_____________________________________________________________________________
331{
333
334 for (G4int i = 0; i < kNofParticlesWSP; i++) {
335 if ((*fIsCutVector)[i]) return true;
336 }
337
338 return false;
339}
340
341//_____________________________________________________________________________
343{
347
348 for (G4int i = 0; i < kNofParticlesWSP; i++) {
349 if ((*fIsControlVector)[i]) return true;
350 }
351
352 return false;
353}
354
355//_____________________________________________________________________________
357{
359
360 for (G4int i = 0; i < kNoG3Controls; i++) {
361 if ((*fControlVector)[i] != kUnsetControlValue) return true;
362 }
363
364 return false;
365}
366
367//_____________________________________________________________________________
369 G4ParticleDefinition* particle) const
370{
373
374 G4String name = particle->GetParticleName();
375 G4String pType = particle->GetParticleType();
376
377 if (name == "gamma") {
378 return kGamma;
379 }
380 else if (name == "e-") {
381 return kElectron;
382 }
383 else if (name == "e+") {
384 return kEplus;
385 }
386 else if ((pType == "baryon" || pType == "meson" || pType == "nucleus" ||
387 pType == "Ion")) {
388 if (particle->GetPDGCharge() == 0) {
389 return kNeutralHadron;
390 }
391 else
392 return kChargedHadron;
393 }
394 else if (name == "mu-" || name == "mu+") {
395 return kMuon;
396 }
397 else {
398 return kNofParticlesWSP;
399 }
400}
401
402//_____________________________________________________________________________
403G4String TG4G3PhysicsManager::GetG3ParticleWSPName(G4int particleWSP) const
404{
407
408 switch (particleWSP) {
409 case kGamma:
410 return "Gamma";
411 break;
412 case kElectron:
413 return "Electron";
414 break;
415 case kEplus:
416 return "Eplus";
417 break;
418 case kNeutralHadron:
419 return "NeutralHadron";
420 break;
421 case kChargedHadron:
422 return "ChargedHadron";
423 break;
424 case kMuon:
425 return "Muon";
426 break;
427 case kAny:
428 return "Any";
429 break;
430 case kNofParticlesWSP:
431 return "NoSP";
432 break;
433 default:
435 "TG4G3PhysicsManager", "GetG3ParticleWSPName", "Wrong particleWSP.");
436 return "";
437 }
438}
Definition of the TG4G3ControlVector class.
Definition of the TG4G3CutVector class.
Definition of the TG4G3Defaults class.
Definition of the TG4G3PhysicsManager class.
Definition of the TG4G3Units class.
Vector of control process values with convenient set/get methods.
static TG4G3Control GetControl(const G4String &controlName)
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.
static TG4G3Cut GetCut(const G4String &cutName)
static G4double Tolerance()
static G4bool CheckCutValue(TG4G3Cut cut, G4double value)
void SetCut(TG4G3Cut cut, G4double cutValue)
Provides a Geant3 way control to Geant4 physics.
G4bool IsGlobalSpecialControls() const
void SwitchIsCutVector(TG4G3Cut cut)
G4bool CheckControlWithTheVector(G4String name, G4double value, TG4G3Control &control, TG4G3ControlValue &controlValue)
TG4G3ControlVector * fControlVector
TG4G3ControlVector.
void SwitchIsControlVector(TG4G3Control control)
G4String GetG3ParticleWSPName(G4int particleWSP) const
void SetCut(TG4G3Cut cut, G4double cutValue)
TG4boolVector * fIsControlVector
vector of booleans which controls are set
TG4boolVector * fIsCutVector
vector of booleans which cuts are set
G4bool fLock
if true: cut/control vectors cannot be modified
TG4G3ParticleWSP GetG3ParticleWSP(G4ParticleDefinition *particle) const
void SetProcess(TG4G3Control control, TG4G3ControlValue controlValue)
static TG4G3PhysicsManager * fgInstance
this instance
TG4G3CutVector * fCutVector
TG4G3CutVector.
G4bool CheckCutWithTheVector(G4String name, G4double value, TG4G3Cut &cut)
static G4double Energy()
Definition TG4G3Units.h:105
static G4double Time()
Definition TG4G3Units.h:93
static void Exception(const TString &className, const TString &methodName, const TString &text)
static TString Endl()
Definition TG4Globals.h:100
TG4G3ControlValue
Enumeration for G3 processes control values.
TG4G3Cut
Enumeration for G3 types of kinetic energy cuts.
Definition TG4G3Cut.h:22
std::vector< G4bool > TG4boolVector
Definition TG4Globals.h:37
TG4G3Control
Enumeration for G3 types of physics processes controls.
@ kUnsetControlValue
value not set
@ kCUTGAM
Definition TG4G3Cut.h:27
@ kDCUTM
Definition TG4G3Cut.h:67
@ kDCUTE
Definition TG4G3Cut.h:62
@ kCUTHAD
Definition TG4G3Cut.h:42
@ kNoG3Cuts
Invalid value.
Definition TG4G3Cut.h:80
@ kBCUTM
Definition TG4G3Cut.h:57
@ kCUTNEU
Definition TG4G3Cut.h:37
@ kCUTMUO
Definition TG4G3Cut.h:47
@ kCUTELE
Definition TG4G3Cut.h:32
@ kPPCUTM
Definition TG4G3Cut.h:72
@ kBCUTE
Definition TG4G3Cut.h:52
@ kTOFMAX
Definition TG4G3Cut.h:77
@ kDCAY
@ kPAIR
@ kPFIS
@ kNoG3Controls
No process control.
@ kHADR
@ kANNI
@ kDRAY
@ kG3LOSS
@ kCOMP
@ kBREM
@ kMUNU
@ kMULS
@ kPHOT
TG4G3ParticleWSP
The particles types which a special process (cuts, controls) is applicable for.
@ kNofParticlesWSP
not a particle with a special control
@ kNeutralHadron
kHADR
@ kGamma
kPAIR, kCOMP, kPHOT, kPHIS
@ kElectron
kDRAY, kBREM, kMULS, kG3LOSS
@ kEplus
kDRAY, kBREM, kMULS, kG3LOSS, kANNI
@ kAny
kDCAY
@ kMuon
kDRAY, kBREM, kMULS, kG3LOSS, kMUNU
@ kChargedHadron
kDRAY, kMULS, kG3LOSS, kHADR,