Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4SDServices.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 "TG4SDServices.h"
16#include "TG4GeometryServices.h"
17#include "TG4Globals.h"
19
20#include <G4LogicalVolume.hh>
21#include <G4LogicalVolumeStore.hh>
22#include <G4Material.hh>
23#include <G4VSensitiveDetector.hh>
24
25#include <TVirtualMCSensitiveDetector.h>
26
27#include <iomanip>
28
31
32G4ThreadLocal std::set<TVirtualMCSensitiveDetector*>* TG4SDServices::fgUserSDs =
33 0;
34G4ThreadLocal std::map<G4String, TVirtualMCSensitiveDetector*>*
36
37//_____________________________________________________________________________
39 : fIsStopRun(false),
40 fVolNameToIdMap(),
41 fVolIdToLVMap(),
42 fLVToVolIdMap(),
43 fIsUserSDs(false)
44{
46
47 if (fgInstance) {
48 TG4Globals::Exception("TG4SDServices", "TG4SDServices",
49 "Cannot create two instances of singleton.");
50 }
51 fgInstance = this;
52}
53
54//_____________________________________________________________________________
60
61//
62// public methods
63//
64
65//_____________________________________________________________________________
67 G4LogicalVolume* lv, G4int id, G4bool fillLVToVolIdMap)
68{
71
72 // cut copy number from sdName
73 G4String volName =
75
76 if (fVolNameToIdMap.find(volName) == fVolNameToIdMap.end())
77 fVolNameToIdMap[volName] = id;
78
79 if (fVolIdToLVMap.find(id) == fVolIdToLVMap.end()) fVolIdToLVMap[id] = lv;
80
81 if (fillLVToVolIdMap) {
82 if (fLVToVolIdMap.find(lv) == fLVToVolIdMap.end()) fLVToVolIdMap[lv] = id;
83 }
84}
85
86//_____________________________________________________________________________
88 const G4String& volumeName, TVirtualMCSensitiveDetector* userSD)
89{
92
93 // G4cout << "TG4SDServices::MapUserSD "
94 // << volumeName << " " << userSD << G4endl;
95
96 // Create the map if it does not yet exist
97 if (!fgUserSDs) {
98 fgUserSDs = new std::set<TVirtualMCSensitiveDetector*>();
99 fgUserSDMap = new std::map<G4String, TVirtualMCSensitiveDetector*>();
100 fIsUserSDs = true;
101 }
102
103 // Add the userSD in the vector only once
104 if (fgUserSDs->find(userSD) == fgUserSDs->end()) {
105 fgUserSDs->insert(userSD);
106 }
107
108 if (fgUserSDMap->find(volumeName) == fgUserSDMap->end()) {
109 (*fgUserSDMap)[volumeName] = userSD;
110 }
111 else {
112 TG4Globals::Warning("TG4SDServices", "MapUserSD",
113 TString("A sensitive detector for volume ") + TString(volumeName.data()) +
114 TString(" has been already defined.") + TG4Globals::Endl() +
115 TString("Setting was ingored."));
116 }
117}
118
119//_____________________________________________________________________________
120void TG4SDServices::PrintStatistics(G4bool open, G4bool close) const
121{
123
124 if (open) TG4Globals::PrintStars(true);
125
126 G4cout << " " << std::setw(5) << NofSensitiveDetectors()
127 << " sensitive detectors" << G4endl;
128
129 if (close) TG4Globals::PrintStars(false);
130}
131
132//_____________________________________________________________________________
134{
136
137 if (fVolNameToIdMap.size()) {
138 G4cout << "Dump of VolNameToIdMap - " << fVolNameToIdMap.size()
139 << " entries:" << G4endl;
140 G4int counter = 0;
141 std::map<G4String, G4int>::const_iterator it;
142 for (it = fVolNameToIdMap.begin(); it != fVolNameToIdMap.end(); ++it) {
143 G4cout << "Map element " << std::right << std::setw(3) << counter++
144 << " ";
145 G4cout << std::left << std::setw(20) << it->first << " ";
146 G4cout << std::right << std::setw(8) << it->second << " ";
147 G4cout << G4endl;
148 }
149 }
150}
151
152//_____________________________________________________________________________
154{
156
157 if (fVolIdToLVMap.size()) {
158 G4cout << "Dump of VolIdToLVMap - " << fVolIdToLVMap.size()
159 << " entries:" << G4endl;
160 G4int counter = 0;
161 std::map<G4int, G4LogicalVolume*>::const_iterator it;
162 for (it = fVolIdToLVMap.begin(); it != fVolIdToLVMap.end(); ++it) {
163 G4cout << "Map element " << std::right << std::setw(3) << counter++
164 << " ";
165 G4cout << std::right << std::setw(8) << it->first << " ";
166 G4cout << std::right << std::setw(12) << it->second << " ";
167 G4cout << std::left << std::setw(20) << it->second->GetName() << " ";
168 G4cout << std::right << G4endl;
169 }
170 }
171}
172
173//_____________________________________________________________________________
175{
177
178 G4LogicalVolumeStore* lvStore = G4LogicalVolumeStore::GetInstance();
179
180 G4cout << "Sensitive volumes (volId, volName): " << G4endl;
181 for (G4int i = 0; i < G4int(lvStore->size()); i++) {
182 G4LogicalVolume* lv = (*lvStore)[i];
183 if (lv->GetSensitiveDetector()) {
184 G4cout << " ";
185 G4cout << std::right << std::setw(4) << GetVolumeID(lv) << " ";
186 G4cout << std::left << std::setw(20) << lv->GetName();
187 G4cout << std::right << G4endl;
188 }
189 }
190}
191
192//_____________________________________________________________________________
194{
196
197 G4LogicalVolumeStore* lvStore = G4LogicalVolumeStore::GetInstance();
198
199 G4cout << "User sensitive detectors (volId, volName, sdName): " << G4endl;
200 for (G4int i = 0; i < G4int(lvStore->size()); i++) {
201 G4LogicalVolume* lv = (*lvStore)[i];
202 TVirtualMCSensitiveDetector* userSD =
203 static_cast<TG4SensitiveDetector*>(lv->GetSensitiveDetector())
204 ->GetUserSD();
205 if (userSD) {
206 G4cout << " ";
207 G4cout << std::right << std::setw(4) << GetVolumeID(lv) << " ";
208 G4cout << std::left << std::setw(20) << lv->GetName() << " ";
209 G4cout << std::left << std::setw(20) << userSD->GetName() << " ";
210 G4cout << std::right << G4endl;
211 }
212 }
213}
214
215//_____________________________________________________________________________
216G4int TG4SDServices::GetVolumeID(const G4String& volName) const
217{
219
220 G4String g4VolName = TG4GeometryServices::Instance()->CutName(volName);
221
222 std::map<G4String, G4int>::const_iterator it = fVolNameToIdMap.find(volName);
223
224 if (it == fVolNameToIdMap.end()) {
225 TG4Globals::Warning("TG4SDServices", "GetVolumeID",
226 "Unknown Volume Id for " + TString(volName.data()));
227 return 0;
228 }
229
230 return it->second;
231}
232
233//_____________________________________________________________________________
234G4int TG4SDServices::GetVolumeID(G4LogicalVolume* logicalVolume) const
235{
239
240 if (fIsUserSDs) {
241 return logicalVolume->GetInstanceID() + fgkFirstVolumeId;
242 }
243
244#ifdef MCDEBUG
245 G4VSensitiveDetector* sd = logicalVolume->GetSensitiveDetector();
246
247 if (sd) {
248 return GetSensitiveDetector(sd)->GetID();
249 }
250 else {
251 std::map<G4LogicalVolume*, G4int>::const_iterator it;
252 it = fLVToVolIdMap.find(logicalVolume);
253 if (it != fLVToVolIdMap.end()) {
254 return it->second;
255 }
256 else {
257 TG4Globals::Warning("TG4SDServices", "GetVolumeID",
258 "Unknown Volume Id for " + TString(logicalVolume->GetName()));
259 return 0;
260 }
261 }
262#else
263 G4VSensitiveDetector* sd = logicalVolume->GetSensitiveDetector();
264 if (sd)
265 return ((TG4SensitiveDetector*)sd)->GetID();
266 else {
267 std::map<G4LogicalVolume*, G4int>::const_iterator it;
268 it = fLVToVolIdMap.find(logicalVolume);
269 if (it != fLVToVolIdMap.end()) {
270 return it->second;
271 }
272 else {
273 TG4Globals::Warning("TG4SDServices", "GetVolumeID",
274 "Unknown Volume Id for" + TString(logicalVolume->GetName()));
275 return 0;
276 }
277 }
278#endif
279}
280
281//_____________________________________________________________________________
282G4int TG4SDServices::GetMediumID(G4LogicalVolume* logicalVolume) const
283{
285
286 if (fIsUserSDs) {
287 return TG4GeometryServices::Instance()->GetMediumId(logicalVolume);
288 }
289
290#ifdef MCDEBUG
291 G4VSensitiveDetector* sd = logicalVolume->GetSensitiveDetector();
292
293 if (sd) {
294 return GetSensitiveDetector(sd)->GetMediumID();
295 }
296 else {
297 return TG4GeometryServices::Instance()->GetMediumId(logicalVolume);
298 }
299#else
300 G4VSensitiveDetector* sd = logicalVolume->GetSensitiveDetector();
301 if (sd)
302 return ((TG4SensitiveDetector*)sd)->GetMediumID();
303 else {
304 return TG4GeometryServices::Instance()->GetMediumId(logicalVolume);
305 }
306#endif
307}
308
309//_____________________________________________________________________________
310TVirtualMCSensitiveDetector* TG4SDServices::GetUserSD(
311 G4String volumeName, G4bool warn) const
312{
313
314 if (!fgUserSDMap) return 0;
315
316 std::map<G4String, TVirtualMCSensitiveDetector*>::const_iterator it =
317 fgUserSDMap->find(volumeName);
318
319 if (it == fgUserSDMap->end()) {
320 if (warn) {
321 TG4Globals::Warning("TG4SDServices", "GetLogicalVolume",
322 "No sensitive detector is defined for volume " +
323 TString(volumeName.data()));
324 }
325 return 0;
326 }
327
328 return it->second;
329}
330
331//_____________________________________________________________________________
332G4String TG4SDServices::GetVolumeName(G4int volumeId) const
333{
337
338 G4LogicalVolume* lv = GetLogicalVolume(volumeId, false);
339
340 if (!lv) {
341 TString text = "volumeId=";
342 text += volumeId;
343 TG4Globals::Warning("TG4SDServices", "GetVolumeName",
344 "Sensitive detector with " + text + " is not defined.");
345 return "";
346 }
347
348 return TG4GeometryServices::Instance()->UserVolumeName(lv->GetName());
349}
350
351//_____________________________________________________________________________
353 G4int volumeId, G4bool warn) const
354{
357
358 std::map<G4int, G4LogicalVolume*>::const_iterator it =
359 fVolIdToLVMap.find(volumeId);
360
361 if (it == fVolIdToLVMap.end()) {
362 if (warn) {
363 TString text = "volumeId=";
364 text += volumeId;
365 TG4Globals::Warning("TG4SDServices", "GetLogicalVolume",
366 "Sensitive detector with " + text + " is not defined.");
367 }
368 return 0;
369 }
370
371 return it->second;
372}
373
374//_____________________________________________________________________________
375Int_t TG4SDServices::GetMediumId(G4int volumeId) const
376{
378
380 GetLogicalVolume(volumeId));
381}
382
383//_____________________________________________________________________________
390
391//_____________________________________________________________________________
393 G4VSensitiveDetector* sd) const
394{
396
397 if (!sd) return 0;
398
399 TG4SensitiveDetector* tsd = dynamic_cast<TG4SensitiveDetector*>(sd);
400
401 if (!tsd) {
402 TG4Globals::Exception("TG4SDServices", "GetSensitiveDetector",
403 "Wrong sensitive detector type.");
404 return 0;
405 }
406 else
407 return tsd;
408}
409
410//_____________________________________________________________________________
411Int_t TG4SDServices::NofVolDaughters(const char* volName) const
412{
414
415 G4int volId = GetVolumeID(volName);
416 G4LogicalVolume* lv = GetLogicalVolume(volId);
417
418 if (!lv) return 0;
419
420 return lv->GetNoDaughters();
421}
422
423//_____________________________________________________________________________
424const char* TG4SDServices::VolDaughterName(const char* volName, Int_t i) const
425{
427
428 G4int volId = GetVolumeID(volName);
429 G4LogicalVolume* lv = GetLogicalVolume(volId);
430
431 if (!lv) return "";
432
433 G4int nofDaughters = lv->GetNoDaughters();
434 if (i < 0 || i >= nofDaughters) {
435 TString text = "index=";
436 text += i;
437 TG4Globals::Warning("TG4SDServices", "VolDaughterName",
438 "Mother volume " + TString(volName) + " has no daughter with " + text +
439 ".");
440 return "";
441 }
442
443 G4VPhysicalVolume* daughter = lv->GetDaughter(i);
444 G4String g4Name = daughter->GetLogicalVolume()->GetName();
445
447}
448
449//_____________________________________________________________________________
450Int_t TG4SDServices::VolDaughterCopyNo(const char* volName, Int_t i) const
451{
453
454 G4int volId = GetVolumeID(volName);
455 G4LogicalVolume* lv = GetLogicalVolume(volId);
456
457 if (!lv) return 0;
458
459 G4int nofDaughters = lv->GetNoDaughters();
460 if (i < 0 || i >= nofDaughters) {
461 TString text = "index=";
462 text += i;
463 TG4Globals::Warning("TG4SDServices", "VolDaughterCopyNo",
464 "Mother volume " + TString(volName) + " has no daughter with " + text +
465 ".");
466 return 0;
467 }
468
469 return lv->GetDaughter(i)->GetCopyNo();
470}
Definition of the TG4GeometryServices class.
Definition of the TG4Globals class and basic container types.
Definition of the TG4SDServices class.
Definition of the TG4SensitiveDetector class.
G4int GetMediumId(G4LogicalVolume *lv) const
const G4String & UserVolumeName(const G4String &name) const
G4String CutName(const char *name) const
static TG4GeometryServices * Instance()
static void PrintStars(G4bool emptyLineFirst)
static void Warning(const TString &className, const TString &methodName, const TString &text)
static void Exception(const TString &className, const TString &methodName, const TString &text)
static TString Endl()
Definition TG4Globals.h:100
Sensitive detectors services.
G4int GetMediumID(G4LogicalVolume *volume) const
static G4ThreadLocal std::set< TVirtualMCSensitiveDetector * > * fgUserSDs
vector of user SDs
static TG4SDServices * fgInstance
this instance
G4bool fIsUserSDs
info about user SDs
static G4ThreadLocal std::map< G4String, TVirtualMCSensitiveDetector * > * fgUserSDMap
map volume name -> user SD
std::map< G4int, G4LogicalVolume * > fVolIdToLVMap
map volume id -> logical volume
TG4SensitiveDetector * GetSensitiveDetector(G4VSensitiveDetector *sd) const
std::map< G4LogicalVolume *, G4int > fLVToVolIdMap
map logical volume -> volume id
TVirtualMCSensitiveDetector * GetUserSD(G4String volumeName, G4bool warn=true) const
void PrintUserSensitiveDetectors() const
G4int GetVolumeID(const G4String &volumeName) const
void PrintSensitiveVolumes() const
void PrintVolNameToIdMap() const
Int_t NofVolDaughters(const char *volName) const
void MapVolume(G4LogicalVolume *lv, G4int id, G4bool fillLVToVolIdMap)
std::map< G4String, G4int > fVolNameToIdMap
map volume name -> volume id
static const G4int fgkFirstVolumeId
the first volume id
Int_t NofSensitiveDetectors() const
void MapUserSD(const G4String &volumeName, TVirtualMCSensitiveDetector *userSD)
void PrintStatistics(G4bool open, G4bool close) const
void PrintVolIdToLVMap() const
Int_t VolDaughterCopyNo(const char *volName, Int_t i) const
G4String GetVolumeName(G4int volumeId) const
G4int GetMediumId(G4int volumeId) const
G4LogicalVolume * GetLogicalVolume(G4int volumeId, G4bool warn=true) const
const char * VolDaughterName(const char *volName, Int_t i) const
Sensitive detector class for calling a user defined stepping function.
static G4int GetTotalNofSensitiveDetectors()
TVirtualMCSensitiveDetector * GetUserSD() const