GGS(GenericGEANT4Simulation)Software  2.99.0
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations 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 232 of file GGSTHitsReader.cpp.

232  {
233  int type = -1;
234  if (_chain) {
235  std::string detectorString{detector};
236  detectorString += "-GGSIntHits";
237  type = _CheckBranch(detectorString.c_str());
238  }
239  if (type == 0 && particleHits) {
240  std::string detectorString{detector};
241  detectorString += "-GGSPartHits";
242  type = _CheckBranch(detectorString.c_str());
243  }
244  if (type == 1 && positionHits) {
245  std::string detectorString{detector};
246  detectorString += "-GGSPosHits";
247  type = _CheckBranch(detectorString.c_str());
248  }
249 
250  return (type != -1);
251 }
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 384 of file GGSTHitsReader.cpp.

384  {
385  for (UInt_t i = 0; i < _detectors.size(); i++) {
386  if (!strcmp(_detectors[i].name.data(), detector))
387  return i;
388  }
389  return -1;
390 }
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 80 of file GGSTHitsReader.cpp.

80  {
81  static const std::string routineName("GGSTHitsReader::GetEntry");
82 
83  if (entry != _chain->GetReadEntry()) {
84  _chain->GetEntry(entry);
85  }
86 
87  if (_chain->GetChainEntryNumber(0) != _firstEntryOfCurrentFile) {
88  GGSCOUT(DEBUG) << "Switching to file " << _chain->GetFile()->GetName() << GGSENDL;
89  _firstEntryOfCurrentFile = _chain->GetChainEntryNumber(0);
90 
91  // Read the GGSHitDetInfo TClonesArray from the new file
92  delete _detInfo;
93  _detInfo = (TClonesArray *)(_chain->GetFile()->Get("GGSHitDetInfo"));
94 
95  // Update the detector indexes in _detectors
96  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++)
97  _detectors[iDet].detInfoArrayIndex = _GetDetectorIndexFromStorage(_detectors[iDet].logVolName.data());
98 
99  // Reload time bins and thresholds
100  _timeBins.clear();
101  _thresholds.clear();
102  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++) {
103  delete _timeBins[iDet];
104  _timeBins.push_back((TArrayF *)(_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-TimeBins")));
105  delete _thresholds[iDet];
106  _thresholds.push_back(
107  (std::vector<float> *)(_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-Thresholds")));
108  }
109  }
110 }
#define GGSENDL
Definition: GGSSmartLog.h:131
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 324 of file GGSTHitsReader.cpp.

324  {
325  int iDet = GetDetectorIndex(detector);
326  if (iDet != -1)
327  return GetHit(iDet, iHit);
328  else
329  return nullptr;
330 }
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 332 of file GGSTHitsReader.cpp.

332  {
333  static TClonesArray *hitArray = nullptr;
334  hitArray = _GetHitArray(iDet);
335  if (hitArray) {
336  GGSTIntHitBase *hit = (GGSTIntHitBase *)(hitArray->At(iHit));
337  if (hit) {
338  // Compute non-persistent informations
339  hit->eDep = hit->eDepTimeBin.GetSum();
340  hit->_hitDetInfo = _GetDetectorInfo(_detectors[iDet].detInfoArrayIndex);
341  hit->_timeBins = _timeBins[iDet];
342  if (hit->_hitDetInfo)
343  hit->_hitVolInfo = (GGSTHitVolInfo *)(hit->_hitDetInfo->volumes.At(hit->_volumeIndex));
344  else
345  hit->_hitVolInfo = nullptr;
346  if (_readPartHits[iDet]) {
347  hit->_partHits = (TClonesArray *)(_partHitArrays[iDet]->At(iHit));
348  }
349  if (_readPosHits[iDet]) {
350  // Set pointer to position hits collection for each particle hit
351  TClonesArray *partHits = hit->_partHits;
352  int nPartHits = partHits->GetEntries();
353  for (int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
354  GGSTPartHit *partHit = (GGSTPartHit *)(partHits->At(iPartHit));
355  partHit->_posHits = (TClonesArray *)(_posHitArrays[iDet]->At(partHit->_posHitIndex));
356  }
357  }
358  return hit;
359  }
360  }
361  return nullptr;
362 }
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 364 of file GGSTHitsReader.cpp.

364  {
365  int iDet = GetDetectorIndex(detector);
366  if (iDet >= 0) {
367  return GetIntHitsClassName(iDet);
368  } else {
369  return "";
370  }
371 }
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 373 of file GGSTHitsReader.cpp.

373  {
374  TClonesArray *hitArray = _GetHitArray(iDet);
375  if (hitArray) {
376  std::string retValue = hitArray->GetName();
377  return retValue.substr(0, retValue.size() - 1);
378  } else {
379  return "";
380  }
381 }
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 422 of file GGSTHitsReader.cpp.

422  {
423  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
424  int detIndex = GetDetectorIndex(detector.c_str());
425  if (detIndex == -1) {
426  throw std::range_error(routineName + ": detector " + detector + " not found");
427  }
428 
429  if (_thresholds[detIndex]) {
430  return (*(_thresholds[detIndex]))[0];
431  } else {
432  throw std::runtime_error(routineName + ": thresholds for integrated hits not available for detector " + detector);
433  }
434 }
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 303 of file GGSTHitsReader.cpp.

303  {
304  static short int iHitArray;
305  iHitArray = GetDetectorIndex(detector);
306  if (iHitArray != -1)
307  return (_GetHitArray(iHitArray))->GetEntries();
308  else
309  return 0;
310 }
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 314 of file GGSTHitsReader.cpp.

314  {
315  TClonesArray *detArray = _GetHitArray(iDet);
316  if (detArray)
317  return detArray->GetEntries();
318  else
319  return 0;
320 }
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 438 of file GGSTHitsReader.cpp.

438  {
439  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
440  int detIndex = GetDetectorIndex(detector.c_str());
441  if (detIndex == -1) {
442  throw std::range_error(routineName + ": detector " + detector + " not found");
443  }
444 
445  if (_thresholds[detIndex]) {
446  return (*(_thresholds[detIndex]))[1];
447  } else {
448  throw std::runtime_error(routineName + ": thresholds for particle hits not available for detector " + detector);
449  }
450 }
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 454 of file GGSTHitsReader.cpp.

454  {
455  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
456  int detIndex = GetDetectorIndex(detector.c_str());
457  if (detIndex == -1) {
458  throw std::range_error(routineName + ": detector " + detector + " not found");
459  }
460 
461  if (_thresholds[detIndex]) {
462  return (*(_thresholds[detIndex]))[2];
463  } else {
464  throw std::runtime_error(routineName + ": thresholds for position hits not available for detector " + detector);
465  }
466 }
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 394 of file GGSTHitsReader.cpp.

394  {
395  if (_chain) {
396  TString branchName(detector);
397  branchName.Append("-GGSPartHits");
398  TBranch *branch = _chain->GetBranch(branchName);
399  if (branch)
400  return true;
401  else
402  return false;
403  } else
404  return false;
405 }
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 409 of file GGSTHitsReader.cpp.

409  {
410  if (_chain) {
411  TString branchName(detector);
412  branchName.Append("-GGSPosHits");
413  TBranch *branch = _chain->GetBranch(branchName);
414  if (branch)
415  return true;
416  else
417  return false;
418  } else
419  return false;
420 }
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  delete _detInfo;
49  _detInfo = (TClonesArray *)(_chain->GetFile()->Get("GGSHitDetInfo"));
50 
51  // Search for hit branches in master chain
52  for (int iBranch = 0; iBranch < _chain->GetListOfBranches()->GetEntries(); iBranch++) {
53  if (TString(((TBranch *)(_chain->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
54  return true;
55  }
56  }
57 
58  // Search for hit branches in friend chains
59  if (_chain->GetListOfFriends()) {
60  for (int iFriend = 0; iFriend < _chain->GetListOfFriends()->GetEntries(); iFriend++) {
61  TTree *currFriend = ((TFriendElement *)(_chain->GetListOfFriends()->At(iFriend)))->GetTree();
62  for (int iBranch = 0; iBranch < currFriend->GetListOfBranches()->GetEntries(); iBranch++) {
63  if (TString(((TBranch *)(currFriend->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
64  return true;
65  }
66  }
67  }
68  }
69 
70  _chain = nullptr;
71  delete _detInfo;
72  _detInfo = nullptr;
73  _firstEntryOfCurrentFile = -1;
74 
75  return false;
76 }
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 114 of file GGSTHitsReader.cpp.

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

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