GGS(GenericGEANT4Simulation)Software  2.6.3
 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)
 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...
 
GGSTIntHitGetHit (const char *detector, unsigned int iHit)
 Get the specified hit. More...
 
GGSTIntHitGetHit (int iDet, unsigned int iHit)
 Get the specified hit. 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)

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.
Returns
true if hits for the specified detector exist in the hits tree.
See Also
SetDetector

Definition at line 237 of file GGSTHitsReader.cpp.

237  {
238 
239  if (_chain) {
240  TString detectorString(detector);
241  detectorString.Append("-GGSIntHits");
242  return _CheckBranch(detectorString.Data());
243  }
244 
245  return false;
246 
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 371 of file GGSTHitsReader.cpp.

371  {
372  for (UInt_t i = 0; i < _detectors.size(); i++) {
373  if (!strcmp(_detectors[i].name.data(), detector))
374  return i;
375  }
376  return -1;
377 }
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  COUT(DEBUG) << "Switching to file " << _chain->GetFile()->GetName() << ENDL;
89  _firstEntryOfCurrentFile = _chain->GetChainEntryNumber(0);
90 
91  // Read the GGSHitDetInfo TClonesArray from the new file
92  _detInfo = (TClonesArray*) (_chain->GetFile()->Get("GGSHitDetInfo"));
93 
94  // Update the detector indexes in _detectors
95  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++)
96  _detectors[iDet].detInfoArrayIndex = _GetDetectorIndexFromStorage(_detectors[iDet].logVolName.data());
97 
98  // Reload time bins and thresholds
99  _timeBins.clear();
100  _thresholds.clear();
101  for (unsigned int iDet = 0; iDet < _detectors.size(); iDet++) {
102  delete _timeBins[iDet];
103  _timeBins.push_back((TArrayF*) (_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-TimeBins")));
104  delete _thresholds[iDet];
105  _thresholds.push_back(
106  (std::vector<float>*) (_chain->GetFile()->Get(TString(_detectors[iDet].name) + "-GGSIntHits-Thresholds")));
107  }
108  }
109 }
#define ENDL
Definition: GGSSmartLog.h:104
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:77
GGSTIntHit * 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 (NULL if iHit is out of bounds).

Definition at line 328 of file GGSTHitsReader.cpp.

328  {
329  int iDet = GetDetectorIndex(detector);
330  if (iDet != -1)
331  return GetHit(iDet, iHit);
332  else
333  return nullptr;
334 }
GGSTIntHit * GetHit(const char *detector, unsigned int iHit)
Get the specified hit.
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.
GGSTIntHit * 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 (NULL if iHit is out of bounds).

Definition at line 336 of file GGSTHitsReader.cpp.

336  {
337  static TClonesArray *hitArray = nullptr;
338  hitArray = _GetHitArray(iDet);
339  if (hitArray) {
340  static GGSTIntHit *hit = nullptr;
341  hit = (GGSTIntHit*) (hitArray->At(iHit));
342  if (hit) {
343  // Compute non-persistent informations
344  hit->eDep = hit->eDepTimeBin.GetSum();
345  hit->_hitDetInfo = _GetDetectorInfo(_detectors[iDet].detInfoArrayIndex);
346  hit->_timeBins = _timeBins[iDet];
347  if (hit->_hitDetInfo)
348  hit->_hitVolInfo = (GGSTHitVolInfo*) (hit->_hitDetInfo->volumes.At(hit->_volumeIndex));
349  else
350  hit->_hitVolInfo = nullptr;
351  if (_readPartHits[iDet]) {
352  hit->_partHits = (TClonesArray*) (_partHitArrays[iDet]->At(iHit));
353  }
354  if (_readPosHits[iDet]) {
355  // Set pointer to position hits collection for each particle hit
356  TClonesArray *partHits = hit->_partHits;
357  int nPartHits = partHits->GetEntries();
358  for (int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
359  GGSTPartHit *partHit = (GGSTPartHit*) (partHits->At(iPartHit));
360  partHit->_posHits = (TClonesArray*) (_posHitArrays[iDet]->At(partHit->_posHitIndex));
361  }
362  }
363  return hit;
364  }
365  }
366  return nullptr;
367 }
Float_t eDep
Deposited energy (transient).
Definition: GGSTHits.h:164
Class to store G4 particle hits.
Definition: GGSTHits.h:82
TClonesArray volumes
Array of GGSTHitVolInfo objects.
TArrayF eDepTimeBin
Energy deposit in each time bin.
Definition: GGSTHits.h:165
Class to store G4 position hits.
Definition: GGSTHits.h:160
GGSTHitVolInfo.h GGSTHitVolInfo class declaration.
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 411 of file GGSTHitsReader.cpp.

411  {
412  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
413  int detIndex = GetDetectorIndex(detector.c_str());
414  if (detIndex == -1) {
415  throw std::range_error(routineName + ": detector " + detector + " not found");
416  }
417 
418  if (_thresholds[detIndex]) {
419  return (*(_thresholds[detIndex]))[0];
420  }
421  else {
422  throw std::runtime_error(routineName + ": thresholds for integrated hits not available for detector " + detector);
423  }
424 }
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 305 of file GGSTHitsReader.cpp.

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

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

428  {
429  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
430  int detIndex = GetDetectorIndex(detector.c_str());
431  if (detIndex == -1) {
432  throw std::range_error(routineName + ": detector " + detector + " not found");
433  }
434 
435  if (_thresholds[detIndex]) {
436  return (*(_thresholds[detIndex]))[1];
437  }
438  else {
439  throw std::runtime_error(routineName + ": thresholds for particle hits not available for detector " + detector);
440  }
441 }
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 445 of file GGSTHitsReader.cpp.

445  {
446  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
447  int detIndex = GetDetectorIndex(detector.c_str());
448  if (detIndex == -1) {
449  throw std::range_error(routineName + ": detector " + detector + " not found");
450  }
451 
452  if (_thresholds[detIndex]) {
453  return (*(_thresholds[detIndex]))[2];
454  }
455  else {
456  throw std::runtime_error(routineName + ": thresholds for position hits not available for detector " + detector);
457  }
458 }
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 381 of file GGSTHitsReader.cpp.

381  {
382  if (_chain) {
383  TString branchName(detector);
384  branchName.Append("-GGSPartHits");
385  TBranch *branch = _chain->GetBranch(branchName);
386  if (branch)
387  return true;
388  else
389  return false;
390  }
391  else
392  return false;
393 }
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 397 of file GGSTHitsReader.cpp.

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

47  {
48 
49  _chain = hitsChain;
50  _detInfo = (TClonesArray*) (_chain->GetFile()->Get("GGSHitDetInfo"));
51 
52  // Search for hit branches in master chain
53  for (int iBranch = 0; iBranch < _chain->GetListOfBranches()->GetEntries(); iBranch++) {
54  if (TString(((TBranch*) (_chain->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
55  return true;
56  }
57  }
58 
59  // Search for hit branches in friend chains
60  if (_chain->GetListOfFriends()) {
61  for (int iFriend = 0; iFriend < _chain->GetListOfFriends()->GetEntries(); iFriend++) {
62  TTree *currFriend = ((TFriendElement*) (_chain->GetListOfFriends()->At(iFriend)))->GetTree();
63  for (int iBranch = 0; iBranch < currFriend->GetListOfBranches()->GetEntries(); iBranch++) {
64  if (TString(((TBranch*) (currFriend->GetListOfBranches()->At(iBranch)))->GetName()).Contains("GGSIntHits")) {
65  return true;
66  }
67  }
68  }
69  }
70 
71  _chain = nullptr;
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 113 of file GGSTHitsReader.cpp.

113  {
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 
128  if (_CheckBranch(detBranchName)) {
129 
130  int detIndex = GetDetectorIndex(detector); // detIndex == -1 if the detector has not been set
131  TBranch *partHitsBranch = _chain->GetBranch(std::string(detector).append("-GGSPartHits").c_str());
132  TBranch *posHitsBranch = _chain->GetBranch(std::string(detector).append("-GGSPosHits").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  posHitsBranch->SetAddress(&(_posHitArrays.back()));
144  _readPosHits.push_back(true);
145  }
146  else {
147  posHitsBranch->SetAddress(&(_posHitArrays[detIndex]));
148  _readPosHits[detIndex] = true;
149  }
150  branches[2] = posHitsBranch;
151  }
152  else {
153  throw std::runtime_error(std::string("PosHits missing for detector ") + detector);
154  }
155  }
156  else {
157  if (detIndex == -1) {
158  _posHitArrays.push_back(nullptr);
159  _readPosHits.push_back(false);
160  }
161  else {
162  if (posHitsBranch) {
163  posHitsBranch->ResetAddress();
164  }
165  _posHitArrays[detIndex] = nullptr; // FIXME: should the element be deleted before this? Is it owned by the file/branch?
166  _readPosHits[detIndex] = false;
167  }
168  }
169 
170  // Particle hits
171  if (readPartHits) {
172  if (partHitsBranch) {
173  infoString.append(" (+PartHits)");
174  if (detIndex == -1) {
175  _partHitArrays.push_back(nullptr);
176  partHitsBranch->SetAddress(&(_partHitArrays.back()));
177  _readPartHits.push_back(true);
178  }
179  else {
180  partHitsBranch->SetAddress(&(_partHitArrays[detIndex]));
181  _readPartHits[detIndex] = true;
182  }
183  branches[1] = partHitsBranch;
184  }
185  else {
186  throw std::runtime_error(std::string("PartHits missing for detector ") + detector);
187  }
188  }
189  else {
190  if (detIndex == -1) {
191  _partHitArrays.push_back(nullptr);
192  _readPartHits.push_back(false);
193  }
194  else {
195  if (partHitsBranch) {
196  partHitsBranch->ResetAddress();
197  }
198  _partHitArrays[detIndex] = nullptr; // FIXME: should the element be deleted before this? Is it owned by the file/branch?
199  _readPartHits[detIndex] = false;
200  }
201  }
202 
203  if (detIndex == -1) {
204  _intHitArrays.push_back(nullptr);
205  _chain->SetBranchAddress(detBranchName, &(_intHitArrays.back()));
206  branches[0] = _chain->GetBranch(detBranchName);
207  _detInfo = (TClonesArray*) (_chain->GetFile()->Get("GGSHitDetInfo"));
208  _timeBins.push_back((TArrayF*) (_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-TimeBins")));
209  _thresholds.push_back(
210  (std::vector<float>*) (_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-Thresholds")));
211 
212  Detector det;
213  det.name = detector;
214  det.logVolName = detector;
215  det.logVolName = det.logVolName.substr(0, det.logVolName.find_first_of('.')); // remove names of scorer and sensitive detector
216  det.detInfoArrayIndex = _GetDetectorIndexFromStorage(det.logVolName.c_str());
217  _detectors.push_back(det);
218  }
219  COUT(INFO) << infoString << ENDL;
220  }
221  else {
222  throw std::runtime_error(std::string("IntHits missing for detector ") + detector);
223  }
224 
225  // Synchronize the branches that have just been set for reading
226  if (_chain->GetReadEntry() != -1) {
227  for (int iBranch = 0; iBranch < 3; ++iBranch) {
228  if (branches[iBranch]) {
229  branches[iBranch]->GetEntry(_chain->GetReadEntry());
230  }
231  }
232  }
233 }
#define ENDL
Definition: GGSSmartLog.h:104
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:77
int GetDetectorIndex(const char *detector)
Retrieves the index of requested detector.

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