GGS(GenericGEANT4Simulation)Software  2.7.0
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
Data Structures | Public Member Functions
GGSTHitsReader Class Reference

Class for reading output of GGSHitsAction. More...

#include <GGSTHitsReader.h>

Inheritance diagram for GGSTHitsReader:
Inheritance graph
[legend]
Collaboration diagram for GGSTHitsReader:
Collaboration graph
[legend]

Public Member Functions

 GGSTHitsReader ()
 Constructor.
 
 ~GGSTHitsReader ()
 Destructor.
 
bool SetChain (TChain *hitsChain)
 Sets the chain. More...
 
void GetEntry (Long64_t entry)
 
void SetDetector (const char *detector, bool readPartHits=false, bool readPosHits=false)
 Enables reading hit informations of the specified detector. More...
 
bool DetectorExists (const char *detector, bool particleHits=false, bool positionHits=false)
 Checks the existence of a given detector within the hits tree. More...
 
void ListDetectors ()
 Prints a list of detectors which are present in hits tree on standard output.
 
std::vector< TString > GetListOfDetectors ()
 Returns the list of detectors which are present in hits tree on standard output.
 
Int_t GetNHits (const char *detector)
 Gets the number of hits for the specified detector. More...
 
Int_t GetNHits (int iDet)
 Gets the number of hits for the specified detector. More...
 
GGSTIntHitBaseGetHit (const char *detector, unsigned int iHit)
 Get the specified hit. More...
 
GGSTIntHitBaseGetHit (int iDet, unsigned int iHit)
 Get the specified hit. More...
 
std::string GetIntHitsClassName (const char *detector)
 Gets the name of the class of the integrated hits for the given detector. More...
 
std::string GetIntHitsClassName (int iDet)
 Gets the name of the class of the integrated hits for the given detector. More...
 
int GetDetectorIndex (const char *detector)
 Retrieves the index of requested detector. More...
 
bool HasPartHits (const char *detector)
 Check if particle hits data is available for the specified detector. More...
 
bool HasPosHits (const char *detector)
 Check if position hits data is available for the specified detector. More...
 
float GetIntHitThreshold (const std::string &detector)
 Retrieves the thresholds used when saving the integrated hits. More...
 
float GetPartHitThreshold (const std::string &detector)
 Retrieves the thresholds used when saving the particle hits. More...
 
float GetPosHitThreshold (const std::string &detector)
 Retrieves the thresholds used when saving the position hits. More...
 
- Public Member Functions inherited from GGSTChainReader
virtual ~GGSTChainReader ()
 Destructor.
 

Detailed Description

Class for reading output of GGSHitsAction.

This class provides methods to read the information about hits in sensitive detectors as saved by GGSHitsAction. Since the hits tree has a branch for every sensitive detector, the user must specify which branches are to be read by calling SetDetector. The hits can be retrieved using GetHit.

Definition at line 37 of file GGSTHitsReader.h.

Member Function Documentation

bool GGSTHitsReader::DetectorExists ( const char *  detector,
bool  particleHits = false,
bool  positionHits = false 
)

Checks the existence of a given detector within the hits tree.

Checks whether the given detector exists inside the hits tree. Note that this routine does not check whether the detector has been set for hit reading or not.

Parameters
detectorThe detector to search for.
particleHitsif true then the existence of particle hits for the given detector will also be checked.
positionHitsif true then the existence of position hits for the given detector will also be checked.
Returns
true if hits for the specified detector with the specified hits exist in the hits tree.
See Also
SetDetector

Definition at line 228 of file GGSTHitsReader.cpp.

228  {
229  int type = -1;
230  if (_chain) {
231  std::string detectorString{detector};
232  detectorString += "-GGSIntHits";
233  type = _CheckBranch(detectorString.c_str());
234  }
235  if (type == 0 && particleHits) {
236  std::string detectorString{detector};
237  detectorString += "-GGSPartHits";
238  type = _CheckBranch(detectorString.c_str());
239  }
240  if (type == 1 && positionHits) {
241  std::string detectorString{detector};
242  detectorString += "-GGSPosHits";
243  type = _CheckBranch(detectorString.c_str());
244  }
245 
246  return (type != -1);
247 }
int GGSTHitsReader::GetDetectorIndex ( const char *  detector)

Retrieves the index of requested detector.

This index can be used to quickly retrieve hits using the GetHit(int, unsigned int) overload of GetHit.

Parameters
detectorThe name of the detector.
Returns
The index of the desired detector.

Definition at line 380 of file GGSTHitsReader.cpp.

380  {
381  for (UInt_t i = 0; i < _detectors.size(); i++) {
382  if (!strcmp(_detectors[i].name.data(), detector))
383  return i;
384  }
385  return -1;
386 }
void GGSTHitsReader::GetEntry ( Long64_t  entry)
virtual

@ brief Reads the specified entry from the hit branches.

Only detectors specified by SetDetector will be read (and only those hits specified in SetDetector call).

Parameters
entryThe desired entry.

Implements GGSTChainReader.

Definition at line 78 of file GGSTHitsReader.cpp.

78  {
79  static const std::string routineName("GGSTHitsReader::GetEntry");
80 
81  if (entry != _chain->GetReadEntry()) {
82  _chain->GetEntry(entry);
83  }
84 
85  if (_chain->GetChainEntryNumber(0) != _firstEntryOfCurrentFile) {
86  COUT(DEBUG) << "Switching to file " << _chain->GetFile()->GetName() << ENDL;
87  _firstEntryOfCurrentFile = _chain->GetChainEntryNumber(0);
88 
89  // Read the GGSHitDetInfo TClonesArray from the new file
90  _detInfo = (TClonesArray *)(_chain->GetFile()->Get("GGSHitDetInfo"));
91 
92  // Update the detector indexes in _detectors
93  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++)
94  _detectors[iDet].detInfoArrayIndex = _GetDetectorIndexFromStorage(_detectors[iDet].logVolName.data());
95 
96  // Reload time bins and thresholds
97  _timeBins.clear();
98  _thresholds.clear();
99  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++) {
100  delete _timeBins[iDet];
101  _timeBins.push_back((TArrayF *)(_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-TimeBins")));
102  delete _thresholds[iDet];
103  _thresholds.push_back(
104  (std::vector<float> *)(_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-Thresholds")));
105  }
106  }
107 }
#define ENDL
Definition: GGSSmartLog.h:105
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:76
GGSTIntHitBase * GGSTHitsReader::GetHit ( const char *  detector,
unsigned int  iHit 
)

Get the specified hit.

This method searches for the desired hit in the specified detector. The detector will be fetched from the detectors internal index by means of a comparison of the name passed as first argument with the stored detector names. This involves a std::string comparison which might hurt performance when done inside the hit loop. If needed, to improve performance use the GetHit(int, unsigned int) overload.

Parameters
detectorThe desired detector.
iHitThe desired hit of the specified detector
Returns
Pointer to the hit object (nullptr if iHit is out of bounds).

Definition at line 320 of file GGSTHitsReader.cpp.

320  {
321  int iDet = GetDetectorIndex(detector);
322  if (iDet != -1)
323  return GetHit(iDet, iHit);
324  else
325  return nullptr;
326 }
GGSTIntHitBase * GetHit(const char *detector, unsigned int iHit)
Get the specified hit.
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
GGSTIntHitBase * GGSTHitsReader::GetHit ( int  iDet,
unsigned int  iHit 
)

Get the specified hit.

The detector index to be passed as first argument can be retrieved by calling GetDetectorIndex.

Parameters
iDetThe detector index.
iHitThe desired hit of the specified detector
Returns
Pointer to the hit object (nullptr if iHit is out of bounds).

Definition at line 328 of file GGSTHitsReader.cpp.

328  {
329  static TClonesArray *hitArray = nullptr;
330  hitArray = _GetHitArray(iDet);
331  if (hitArray) {
332  GGSTIntHitBase *hit = (GGSTIntHitBase *)(hitArray->At(iHit));
333  if (hit) {
334  // Compute non-persistent informations
335  hit->eDep = hit->eDepTimeBin.GetSum();
336  hit->_hitDetInfo = _GetDetectorInfo(_detectors[iDet].detInfoArrayIndex);
337  hit->_timeBins = _timeBins[iDet];
338  if (hit->_hitDetInfo)
339  hit->_hitVolInfo = (GGSTHitVolInfo *)(hit->_hitDetInfo->volumes.At(hit->_volumeIndex));
340  else
341  hit->_hitVolInfo = nullptr;
342  if (_readPartHits[iDet]) {
343  hit->_partHits = (TClonesArray *)(_partHitArrays[iDet]->At(iHit));
344  }
345  if (_readPosHits[iDet]) {
346  // Set pointer to position hits collection for each particle hit
347  TClonesArray *partHits = hit->_partHits;
348  int nPartHits = partHits->GetEntries();
349  for (int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
350  GGSTPartHit *partHit = (GGSTPartHit *)(partHits->At(iPartHit));
351  partHit->_posHits = (TClonesArray *)(_posHitArrays[iDet]->At(partHit->_posHitIndex));
352  }
353  }
354  return hit;
355  }
356  }
357  return nullptr;
358 }
Class to store G4 particle hits.
Definition: GGSTHits.h:176
TClonesArray volumes
Array of GGSTHitVolInfo objects.
Float_t eDep
Deposited energy (transient).
Definition: GGSTHits.h:222
TArrayF eDepTimeBin
Energy deposit in each time bin.
Definition: GGSTHits.h:223
GGSTHitVolInfo.h GGSTHitVolInfo class declaration.
Base class for storing G4 integrated hits.
Definition: GGSTHits.h:218
std::string GGSTHitsReader::GetIntHitsClassName ( const char *  detector)

Gets the name of the class of the integrated hits for the given detector.

Parameters
detectorThe name of the detector;
Returns
The position hits class name.

Definition at line 360 of file GGSTHitsReader.cpp.

360  {
361  int iDet = GetDetectorIndex(detector);
362  if (iDet >= 0) {
363  return GetIntHitsClassName(iDet);
364  } else {
365  return "";
366  }
367 }
std::string GetIntHitsClassName(const char *detector)
Gets the name of the class of the integrated hits for the given detector.
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
std::string GGSTHitsReader::GetIntHitsClassName ( int  iDet)

Gets the name of the class of the integrated hits for the given detector.

Parameters
iDetThe index of the detector;
Returns
The position hits class name.

Definition at line 369 of file GGSTHitsReader.cpp.

369  {
370  TClonesArray *hitArray = _GetHitArray(iDet);
371  if (hitArray) {
372  std::string retValue = hitArray->GetName();
373  return retValue.substr(0, retValue.size() - 1);
374  } else {
375  return "";
376  }
377 }
float GGSTHitsReader::GetIntHitThreshold ( const std::string &  detector)

Retrieves the thresholds used when saving the integrated hits.

Parameters
detectorThe name of the detector.
Returns
The value of the threshold [GeV].
Exceptions
std::range_errorif the given detector is not present.
std::runtime_errorif thresholds are not available in the input file.

Definition at line 418 of file GGSTHitsReader.cpp.

418  {
419  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
420  int detIndex = GetDetectorIndex(detector.c_str());
421  if (detIndex == -1) {
422  throw std::range_error(routineName + ": detector " + detector + " not found");
423  }
424 
425  if (_thresholds[detIndex]) {
426  return (*(_thresholds[detIndex]))[0];
427  } else {
428  throw std::runtime_error(routineName + ": thresholds for integrated hits not available for detector " + detector);
429  }
430 }
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
Int_t GGSTHitsReader::GetNHits ( const char *  detector)

Gets the number of hits for the specified detector.

Parameters
detectorThe desired detector.
Returns
The number of hits (0 if the specified detector does not exist.

Definition at line 299 of file GGSTHitsReader.cpp.

299  {
300  static short int iHitArray;
301  iHitArray = GetDetectorIndex(detector);
302  if (iHitArray != -1)
303  return (_GetHitArray(iHitArray))->GetEntries();
304  else
305  return 0;
306 }
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
Int_t GGSTHitsReader::GetNHits ( int  iDet)

Gets the number of hits for the specified detector.

The detector index to be passed as first parameter can be retrieved by calling GetDetectorIndex.

Parameters
iDetIndex of the desired detector.
Returns
The number of hits (0 if the specified detector does not exist.

Definition at line 310 of file GGSTHitsReader.cpp.

310  {
311  TClonesArray *detArray = _GetHitArray(iDet);
312  if (detArray)
313  return detArray->GetEntries();
314  else
315  return 0;
316 }
float GGSTHitsReader::GetPartHitThreshold ( const std::string &  detector)

Retrieves the thresholds used when saving the particle hits.

Parameters
detectorThe name of the detector.
Returns
The value of the threshold [GeV].
Exceptions
std::range_errorif the given detector is not present.
std::runtime_errorif thresholds are not available in the input file.

Definition at line 434 of file GGSTHitsReader.cpp.

434  {
435  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
436  int detIndex = GetDetectorIndex(detector.c_str());
437  if (detIndex == -1) {
438  throw std::range_error(routineName + ": detector " + detector + " not found");
439  }
440 
441  if (_thresholds[detIndex]) {
442  return (*(_thresholds[detIndex]))[1];
443  } else {
444  throw std::runtime_error(routineName + ": thresholds for particle hits not available for detector " + detector);
445  }
446 }
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
float GGSTHitsReader::GetPosHitThreshold ( const std::string &  detector)

Retrieves the thresholds used when saving the position hits.

Parameters
detectorThe name of the detector.
Returns
The value of the threshold [GeV].
Exceptions
std::range_errorif the given detector is not present.
std::runtime_errorif thresholds are not available in the input file.

Definition at line 450 of file GGSTHitsReader.cpp.

450  {
451  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
452  int detIndex = GetDetectorIndex(detector.c_str());
453  if (detIndex == -1) {
454  throw std::range_error(routineName + ": detector " + detector + " not found");
455  }
456 
457  if (_thresholds[detIndex]) {
458  return (*(_thresholds[detIndex]))[2];
459  } else {
460  throw std::runtime_error(routineName + ": thresholds for position hits not available for detector " + detector);
461  }
462 }
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
bool GGSTHitsReader::HasPartHits ( const char *  detector)

Check if particle hits data is available for the specified detector.

Parameters
detectorThe name of the detector.
Returns
true if particle hits data is available, false otherwise.

Definition at line 390 of file GGSTHitsReader.cpp.

390  {
391  if (_chain) {
392  TString branchName(detector);
393  branchName.Append("-GGSPartHits");
394  TBranch *branch = _chain->GetBranch(branchName);
395  if (branch)
396  return true;
397  else
398  return false;
399  } else
400  return false;
401 }
bool GGSTHitsReader::HasPosHits ( const char *  detector)

Check if position hits data is available for the specified detector.

Parameters
detectorThe name of the detector.
Returns
true if position hits data is available, false otherwise.

Definition at line 405 of file GGSTHitsReader.cpp.

405  {
406  if (_chain) {
407  TString branchName(detector);
408  branchName.Append("-GGSPosHits");
409  TBranch *branch = _chain->GetBranch(branchName);
410  if (branch)
411  return true;
412  else
413  return false;
414  } else
415  return false;
416 }
bool GGSTHitsReader::SetChain ( TChain *  hitsChain)
virtual

Sets the chain.

This method sets the hits chain, from which the informations will be read.

Parameters
hitsChainThe chain containing the hits.

Implements GGSTChainReader.

Definition at line 46 of file GGSTHitsReader.cpp.

46  {
47  _chain = hitsChain;
48  _detInfo = (TClonesArray *)(_chain->GetFile()->Get("GGSHitDetInfo"));
49 
50  // Search for hit branches in master chain
51  for (int iBranch = 0; iBranch < _chain->GetListOfBranches()->GetEntries(); iBranch++) {
52  if (TString(((TBranch *)(_chain->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
53  return true;
54  }
55  }
56 
57  // Search for hit branches in friend chains
58  if (_chain->GetListOfFriends()) {
59  for (int iFriend = 0; iFriend < _chain->GetListOfFriends()->GetEntries(); iFriend++) {
60  TTree *currFriend = ((TFriendElement *)(_chain->GetListOfFriends()->At(iFriend)))->GetTree();
61  for (int iBranch = 0; iBranch < currFriend->GetListOfBranches()->GetEntries(); iBranch++) {
62  if (TString(((TBranch *)(currFriend->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
63  return true;
64  }
65  }
66  }
67  }
68 
69  _chain = nullptr;
70  _detInfo = nullptr;
71  _firstEntryOfCurrentFile = -1;
72 
73  return false;
74 }
void GGSTHitsReader::SetDetector ( const char *  detector,
bool  readPartHits = false,
bool  readPosHits = false 
)

Enables reading hit informations of the specified detector.

This routine will enable the reading of hits for the specified detector. Hits of other detector will not be read to improve I/O performance. The flags allow for reading particle and/or position hits if available (by default they will not be read from ROOT file, again to maximize I/O performance). Reading of position hits requires reading particle hits, so a true value for readPosHits will force also readPartHits to true. This method throws a std::runtime_exception if the given detector does not exist or if the requested particle and/or position hits are not available, and provides the strong exception guarantee (i.e. no modification to the reader is done if an exception is thrown).

Parameters
detectorThe detector name.
readPartHitsif true, particle hits will be read (if available).
readPosHitsif true, position hits will be read (if available).
Exceptions
std::runtime_errorif the requested detector+hits combination is not available.

Definition at line 111 of file GGSTHitsReader.cpp.

111  {
112  static const std::string routineName("GGSTHitsReader::SetDetector");
113 
114  if (!_chain)
115  return;
116 
117  if (readPosHits) {
118  readPartHits = true;
119  }
120 
121  TString detBranchName = detector;
122  detBranchName.Append("-GGSIntHits");
123  TBranch *branches[3] = {nullptr, nullptr, nullptr}; // New branches (for integrated, particle and position hits)
124  if (_CheckBranch(detBranchName) == 0) {
125  int detIndex = GetDetectorIndex(detector); // detIndex == -1 if the detector has not been set
126  std::string partHitsBranchName{std::string(detector).append("-GGSPartHits")};
127  std::string posHitsBranchName{std::string(detector).append("-GGSPartHits")};
128  TBranch *partHitsBranch = _chain->GetBranch(partHitsBranchName.c_str());
129  TBranch *posHitsBranch = _chain->GetBranch(posHitsBranchName.c_str());
130 
131  std::string infoString("Set detector ");
132  infoString.append(detector);
133 
134  // Position hits
135  if (readPosHits) {
136  if (posHitsBranch) {
137  infoString.append(" (+PosHits)");
138  if (detIndex == -1) {
139  _posHitArrays.push_back(nullptr);
140  _chain->SetBranchAddress(posHitsBranchName.c_str(), &(_posHitArrays.back()));
141  _readPosHits.push_back(true);
142  } else {
143  _chain->SetBranchAddress(posHitsBranchName.c_str(), &(_posHitArrays[detIndex]));
144  _readPosHits[detIndex] = true;
145  }
146  branches[2] = posHitsBranch;
147  } else {
148  throw std::runtime_error(std::string("PosHits missing for detector ") + detector);
149  }
150  } else {
151  if (detIndex == -1) {
152  _posHitArrays.push_back(nullptr);
153  _readPosHits.push_back(false);
154  } else {
155  if (posHitsBranch) {
156  posHitsBranch->ResetAddress();
157  }
158  _posHitArrays[detIndex] = nullptr; // FIXME: should the element be deleted before this? Is it
159  // owned by the file/branch?
160  _readPosHits[detIndex] = false;
161  }
162  }
163 
164  // Particle hits
165  if (readPartHits) {
166  if (partHitsBranch) {
167  infoString.append(" (+PartHits)");
168  if (detIndex == -1) {
169  _partHitArrays.push_back(nullptr);
170  _chain->SetBranchAddress(partHitsBranchName.c_str(), &(_partHitArrays.back()));
171  _readPartHits.push_back(true);
172  } else {
173  _chain->SetBranchAddress(partHitsBranchName.c_str(), &(_partHitArrays[detIndex]));
174  _readPartHits[detIndex] = true;
175  }
176  branches[1] = partHitsBranch;
177  } else {
178  throw std::runtime_error(std::string("PartHits missing for detector ") + detector);
179  }
180  } else {
181  if (detIndex == -1) {
182  _partHitArrays.push_back(nullptr);
183  _readPartHits.push_back(false);
184  } else {
185  if (partHitsBranch) {
186  partHitsBranch->ResetAddress();
187  }
188  _partHitArrays[detIndex] = nullptr; // FIXME: should the element be deleted before this? Is it
189  // owned by the file/branch?
190  _readPartHits[detIndex] = false;
191  }
192  }
193 
194  if (detIndex == -1) {
195  _intHitArrays.push_back(nullptr);
196  _chain->SetBranchAddress(detBranchName, &(_intHitArrays.back()));
197  branches[0] = _chain->GetBranch(detBranchName);
198  _detInfo = (TClonesArray *)(_chain->GetFile()->Get("GGSHitDetInfo"));
199  _timeBins.push_back((TArrayF *)(_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-TimeBins")));
200  _thresholds.push_back(
201  (std::vector<float> *)(_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-Thresholds")));
202 
203  Detector det;
204  det.name = detector;
205  det.logVolName = detector;
206  det.logVolName =
207  det.logVolName.substr(0, det.logVolName.find_first_of('.')); // remove names of scorer and sensitive detector
208  det.detInfoArrayIndex = _GetDetectorIndexFromStorage(det.logVolName.c_str());
209  _detectors.push_back(det);
210  }
211  COUT(INFO) << infoString << ENDL;
212  } else {
213  throw std::runtime_error(std::string("IntHits missing for detector ") + detector);
214  }
215 
216  // Synchronize the branches that have just been set for reading
217  if (_chain->GetReadEntry() != -1) {
218  for (int iBranch = 0; iBranch < 3; ++iBranch) {
219  if (branches[iBranch]) {
220  branches[iBranch]->GetEntry(_chain->GetReadEntry());
221  }
222  }
223  }
224 }
#define ENDL
Definition: GGSSmartLog.h:105
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:76
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.

The documentation for this class was generated from the following files: