GGS(GenericGEANT4Simulation)Software  2.6.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)
 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 207 of file GGSTHitsReader.cpp.

207  {
208 
209  if (_chain) {
210  TString detectorString(detector);
211  detectorString.Append("-GGSIntHits");
212  return _CheckBranch(detectorString.Data());
213  }
214 
215  return false;
216 
217 }
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 338 of file GGSTHitsReader.cpp.

338  {
339  for (UInt_t i = 0; i < _detectors.size(); i++) {
340  if (!strcmp(_detectors[i].name.data(), detector))
341  return i;
342  }
343  return -1;
344 }
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:93
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:66
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 295 of file GGSTHitsReader.cpp.

295  {
296  int iDet = GetDetectorIndex(detector);
297  if (iDet != -1)
298  return GetHit(iDet, iHit);
299  else
300  return NULL;
301 }
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 303 of file GGSTHitsReader.cpp.

303  {
304  static TClonesArray* hitArray = NULL;
305  hitArray = _GetHitArray(iDet);
306  if (hitArray) {
307  static GGSTIntHit *hit = NULL;
308  hit = (GGSTIntHit*) (hitArray->At(iHit));
309  if (hit) {
310  // Compute non-persistent informations
311  hit->eDep = hit->eDepTimeBin.GetSum();
312  hit->_hitDetInfo = _GetDetectorInfo(_detectors[iDet].detInfoArrayIndex);
313  hit->_timeBins = _timeBins[iDet];
314  if (hit->_hitDetInfo)
315  hit->_hitVolInfo = (GGSTHitVolInfo*) (hit->_hitDetInfo->volumes.At(hit->_volumeIndex));
316  else
317  hit->_hitVolInfo = NULL;
318  if (_readPartHits[iDet]) {
319  hit->_partHits = (TClonesArray*) (_partHitArrays[iDet]->At(iHit));
320  }
321  if (_readPosHits[iDet]) {
322  // Set pointer to position hits collection for each particle hit
323  TClonesArray *partHits = hit->_partHits;
324  int nPartHits = partHits->GetEntries();
325  for (int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
326  GGSTPartHit *partHit = (GGSTPartHit*) (partHits->At(iPartHit));
327  partHit->_posHits = (TClonesArray*) (_posHitArrays[iDet]->At(partHit->_posHitIndex));
328  }
329  }
330  return hit;
331  }
332  }
333  return NULL;
334 }
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 378 of file GGSTHitsReader.cpp.

378  {
379  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
380  int detIndex = GetDetectorIndex(detector.c_str());
381  if (detIndex == -1) {
382  throw std::range_error(routineName + ": detector " + detector + " not found");
383  }
384 
385  if (_thresholds[detIndex]) {
386  return (*(_thresholds[detIndex]))[0];
387  }
388  else {
389  throw std::runtime_error(routineName + ": thresholds for integrated hits not available for detector " + detector);
390  }
391 }
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 272 of file GGSTHitsReader.cpp.

272  {
273  static short int iHitArray;
274  iHitArray = GetDetectorIndex(detector);
275  if (iHitArray != -1)
276  return (_GetHitArray(iHitArray))->GetEntries();
277  else
278  return 0;
279 }
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 283 of file GGSTHitsReader.cpp.

283  {
284 
285  TClonesArray *detArray = _GetHitArray(iDet);
286  if (detArray)
287  return detArray->GetEntries();
288  else
289  return 0;
290 
291 }
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 395 of file GGSTHitsReader.cpp.

395  {
396  const std::string routineName("GSTHitsReader::GetIntHitThreshold");
397  int detIndex = GetDetectorIndex(detector.c_str());
398  if (detIndex == -1) {
399  throw std::range_error(routineName + ": detector " + detector + " not found");
400  }
401 
402  if (_thresholds[detIndex]) {
403  return (*(_thresholds[detIndex]))[1];
404  }
405  else {
406  throw std::runtime_error(routineName + ": thresholds for particle hits not available for detector " + detector);
407  }
408 }
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 412 of file GGSTHitsReader.cpp.

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

348  {
349  if (_chain) {
350  TString branchName(detector);
351  branchName.Append("-GGSPartHits");
352  TBranch *branch = _chain->GetBranch(branchName);
353  if (branch)
354  return true;
355  else
356  return false;
357  }
358  else
359  return false;
360 }
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 364 of file GGSTHitsReader.cpp.

364  {
365  if (_chain) {
366  TString branchName(detector);
367  branchName.Append("-GGSPosHits");
368  TBranch *branch = _chain->GetBranch(branchName);
369  if (branch)
370  return true;
371  else
372  return false;
373  }
374  else
375  return false;
376 }
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 = NULL;
72  _detInfo = NULL;
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  if (_CheckDetector(detector) && _CheckBranch(detBranchName)) {
128 
129  std::string infoString("Set detector ");
130  infoString.append(detector);
131 
132  if (readPartHits) {
133  TString branchName(detector);
134  branchName.Append("-GGSPartHits");
135  TBranch *branch = _chain->GetBranch(branchName);
136  if (branch) {
137  infoString.append(" (+PartHits)");
138  _partHitArrays.push_back(NULL);
139  _chain->SetBranchAddress(branchName, &(_partHitArrays.back()));
140  branches[1] = _chain->GetBranch(branchName);
141  _readPartHits.push_back(true);
142  }
143  else {
144  throw std::runtime_error(std::string("PartHits missing for detector ") + detector);
145  }
146  }
147  else {
148  _partHitArrays.push_back(NULL);
149  _readPartHits.push_back(false);
150  }
151 
152  if (readPosHits) {
153  TString branchName(detector);
154  branchName.Append("-GGSPosHits");
155  TBranch *branch = _chain->GetBranch(branchName);
156  if (branch) {
157  infoString.append(" (+PosHits)");
158  _posHitArrays.push_back(NULL);
159  _chain->SetBranchAddress(branchName, &(_posHitArrays.back()));
160  branches[2] = _chain->GetBranch(branchName);
161  _readPosHits.push_back(true);
162  }
163  else {
164  _partHitArrays.pop_back();
165  _readPartHits.pop_back();
166  _chain->ResetBranchAddress(_chain->GetBranch((std::string(detector) + "-GGSPartHits").c_str()));
167  throw std::runtime_error(std::string("PosHits missing for detector ") + detector);
168  }
169  }
170  else {
171  _posHitArrays.push_back(NULL);
172  _readPosHits.push_back(false);
173  }
174 
175  _intHitArrays.push_back(NULL);
176  _chain->SetBranchAddress(detBranchName, &(_intHitArrays.back()));
177  branches[0] = _chain->GetBranch(detBranchName);
178  _detInfo = (TClonesArray*) (_chain->GetFile()->Get("GGSHitDetInfo"));
179  _timeBins.push_back((TArrayF*) (_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-TimeBins")));
180  _thresholds.push_back((std::vector<float>*) (_chain->GetFile()->Get(TString(detector) + "-GGSIntHits-Thresholds")));
181 
182  Detector det;
183  det.name = detector;
184  det.logVolName = detector;
185  det.logVolName = det.logVolName.substr(0, det.logVolName.find_first_of('.')); // remove names of scorer and sensitive detector
186  det.detInfoArrayIndex = _GetDetectorIndexFromStorage(det.logVolName.c_str());
187  _detectors.push_back(det);
188 
189  COUT(INFO) << infoString << ENDL;
190  }
191  else {
192  throw std::runtime_error(std::string("IntHits missing for detector ") + detector);
193  }
194 
195  // Synchronize the branches that have just been set for reading
196  if (_chain->GetReadEntry() != -1) {
197  for (int iBranch = 0; iBranch < 3; ++iBranch) {
198  if (branches[iBranch]) {
199  branches[iBranch]->GetEntry(_chain->GetReadEntry());
200  }
201  }
202  }
203 }
#define ENDL
Definition: GGSSmartLog.h:93
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:66

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