StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StPicoDstReader.cxx
1 //
2 // StPicoDstReader allows to read picoDst file or a list of files
3 //
4 
5 // C++ headers
6 #include <string>
7 #include <sstream>
8 #include <iostream>
9 #include <fstream>
10 #include <assert.h>
11 
12 // PicoDst headers
13 #include "StPicoMessMgr.h"
14 #include "StPicoDstReader.h"
15 #include "StPicoEvent.h"
16 #include "StPicoTrack.h"
17 #include "StPicoEmcTrigger.h"
18 #include "StPicoMtdTrigger.h"
19 #include "StPicoBbcHit.h"
20 #include "StPicoEpdHit.h"
21 #include "StPicoBTowHit.h"
22 #include "StPicoBTofHit.h"
23 #include "StPicoMtdHit.h"
24 #include "StPicoFmsHit.h"
25 #include "StPicoBEmcPidTraits.h"
26 #include "StPicoBTofPidTraits.h"
27 #include "StPicoMtdPidTraits.h"
28 #include "StPicoTrackCovMatrix.h"
29 #include "StPicoBEmcSmdEHit.h"
30 #include "StPicoBEmcSmdPHit.h"
31 #include "StPicoETofHit.h"
32 #include "StPicoETofPidTraits.h"
33 #include "StPicoMcVertex.h"
34 #include "StPicoMcTrack.h"
35 #include "StPicoArrays.h"
36 #include "StPicoDst.h"
37 
38 // ROOT headers
39 #include "TRegexp.h"
40 
41 ClassImp(StPicoDstReader)
42 
43 //_________________
44 StPicoDstReader::StPicoDstReader(const Char_t* inFileName) :
45  mPicoDst(new StPicoDst()), mChain(NULL), mTree(NULL),
46  mEventCounter(0), mPicoArrays{}, mStatusArrays{} {
47 
48  streamerOff();
49  createArrays();
50  std::fill_n(mStatusArrays, sizeof(mStatusArrays) / sizeof(mStatusArrays[0]), 1);
51  mInputFileName = inFileName;
52 }
53 
54 //_________________
56  if(mChain) {
57  delete mChain;
58  }
59  if(mPicoDst) {
60  delete mPicoDst;
61  }
62 }
63 
64 //_________________
65 void StPicoDstReader::clearArrays() {
66  for(Int_t iArr=0; iArr<StPicoArrays::NAllPicoArrays; iArr++) {
67  mPicoArrays[iArr]->Clear();
68  }
69 }
70 
71 //_________________
72 void StPicoDstReader::SetStatus(const Char_t *branchNameRegex, Int_t enable) {
73  if(strncmp(branchNameRegex, "St", 2) == 0) {
74  // Ignore first "St"
75  branchNameRegex += 2;
76  }
77 
78  TRegexp re(branchNameRegex, 1);
79  for(Int_t iArr=0; iArr<StPicoArrays::NAllPicoArrays; iArr++) {
80  Ssiz_t len;
81  if(re.Index(StPicoArrays::picoArrayNames[iArr], &len) < 0) continue;
82  LOG_INFO << "StPicoDstMaker::SetStatus " << enable
83  << " to " << StPicoArrays::picoArrayNames[iArr] << endm;
84  mStatusArrays[iArr] = enable;
85  }
86 
87  setBranchAddresses(mChain);
88 }
89 
90 //_________________
91 void StPicoDstReader::setBranchAddresses(TChain *chain) {
92  if (!chain) return;
93  chain->SetBranchStatus("*", 0);
94  TString ts;
95  for (Int_t i = 0; i < StPicoArrays::NAllPicoArrays; ++i) {
96  if (mStatusArrays[i] == 0) continue;
97  char const* bname = StPicoArrays::picoArrayNames[i];
98  TBranch* tb = chain->GetBranch(bname);
99  if (!tb) {
100  LOG_WARN << "setBranchAddress: Branch name " << bname << " does not exist!" << endm;
101  continue;
102  }
103  ts = bname;
104  ts += "*";
105  chain->SetBranchStatus(ts, 1);
106  chain->SetBranchAddress(bname, mPicoArrays + i);
107  assert(tb->GetAddress() == (char*)(mPicoArrays + i));
108  }
109  mTree = mChain->GetTree();
110 }
111 
112 //_________________
113 void StPicoDstReader::streamerOff() {
114  // This is to to save space on the file. No need for TObject bits for this structure.
115  // see: https://root.cern.ch/doc/master/classTClass.html#a606b0442d6fec4b1cd52f43bca73aa51
116  StPicoEvent::Class()->IgnoreTObjectStreamer();
117  StPicoTrack::Class()->IgnoreTObjectStreamer();
118  StPicoBTofHit::Class()->IgnoreTObjectStreamer();
119  StPicoBTowHit::Class()->IgnoreTObjectStreamer();
120  StPicoMtdHit::Class()->IgnoreTObjectStreamer();
121  StPicoBbcHit::Class()->IgnoreTObjectStreamer();
122  StPicoEpdHit::Class()->IgnoreTObjectStreamer();
123  StPicoFmsHit::Class()->IgnoreTObjectStreamer();
124  StPicoEmcTrigger::Class()->IgnoreTObjectStreamer();
125  StPicoMtdTrigger::Class()->IgnoreTObjectStreamer();
126  StPicoBTofPidTraits::Class()->IgnoreTObjectStreamer();
127  StPicoBEmcPidTraits::Class()->IgnoreTObjectStreamer();
128  StPicoMtdPidTraits::Class()->IgnoreTObjectStreamer();
129  StPicoTrackCovMatrix::Class()->IgnoreTObjectStreamer();
130  StPicoBEmcSmdEHit::Class()->IgnoreTObjectStreamer();
131  StPicoBEmcSmdPHit::Class()->IgnoreTObjectStreamer();
132  StPicoETofHit::Class()->IgnoreTObjectStreamer();
133  StPicoETofPidTraits::Class()->IgnoreTObjectStreamer();
134  StPicoMcVertex::Class()->IgnoreTObjectStreamer();
135  StPicoMcTrack::Class()->IgnoreTObjectStreamer();
136 }
137 
138 //_________________
139 void StPicoDstReader::createArrays() {
140 
141  for(Int_t iArr=0; iArr<StPicoArrays::NAllPicoArrays; iArr++) {
142  mPicoArrays[iArr] = new TClonesArray(StPicoArrays::picoArrayTypes[iArr],
144  }
145  mPicoDst->set(mPicoArrays);
146 }
147 
148 //_________________
150  if(mChain) {
151  delete mChain;
152  }
153  mChain = NULL;
154 }
155 
156 //_________________
158 
159  if(!mChain) {
160  mChain = new TChain("PicoDst");
161  }
162 
163  std::string const dirFile = mInputFileName.Data();
164 
165  if( dirFile.find(".list") != std::string::npos ||
166  dirFile.find(".lis") != std::string::npos ) {
167 
168  std::ifstream inputStream( dirFile.c_str() );
169 
170  if(!inputStream) {
171  LOG_ERROR << "ERROR: Cannot open list file " << dirFile << endm;
172  }
173 
174  Int_t nFile = 0;
175  std::string file;
176  size_t pos;
177  while(getline(inputStream, file)) {
178  // NOTE: our external formatters may pass "file NumEvents"
179  // Take only the first part
180  //cout << "DEBUG found " << file << endl;
181  pos = file.find_first_of(" ");
182  if (pos != std::string::npos ) file.erase(pos,file.length()-pos);
183  //cout << "DEBUG found [" << file << "]" << endl;
184 
185  if(file.find(".picoDst.root") != std::string::npos) {
186  TFile* ftmp = TFile::Open(file.c_str());
187  if(ftmp && !ftmp->IsZombie() && ftmp->GetNkeys()) {
188  LOG_INFO << " Read in picoDst file " << file << endm;
189  mChain->Add(file.c_str());
190  ++nFile;
191  } //if(ftmp && !ftmp->IsZombie() && ftmp->GetNkeys())
192 
193  if (ftmp) {
194  ftmp->Close();
195  } //if (ftmp)
196  } //if(file.find(".picoDst.root") != std::string::npos)
197  } //while (getline(inputStream, file))
198 
199  LOG_INFO << " Total " << nFile << " files have been read in. " << endm;
200  } //if(dirFile.find(".list") != std::string::npos || dirFile.find(".lis" != string::npos))
201  else if(dirFile.find(".picoDst.root") != std::string::npos) {
202  mChain->Add(dirFile.c_str());
203  }
204  else {
205  LOG_WARN << " No good input file to read ... " << endm;
206  }
207 
208  if(mChain) {
209  setBranchAddresses(mChain);
210  mChain->SetCacheSize(50e6);
211  mChain->AddBranchToCache("*");
212  mPicoDst->set(mPicoArrays);
213  }
214 }
215 
216 //_________________
217 Bool_t StPicoDstReader::readPicoEvent(Long64_t iEvent __attribute__((unused)) ) {
218 
219  Int_t mStatusRead = true; // true - okay, false - nothing to read
220 
221  if (!mChain) {
222  LOG_WARN << " No input files ... ! EXIT" << endm;
223  mStatusRead = false;
224  return mStatusRead;
225  }
226 
227  Int_t bytes = mChain->GetEntry(mEventCounter++);
228  Int_t nCycles = 0;
229  while( bytes <= 0) {
230  if( mEventCounter >= mChain->GetEntriesFast() ) {
231  }
232 
233  LOG_WARN << "Encountered invalid entry or I/O error while reading event "
234  << mEventCounter << " from \"" << mChain->GetName() << "\" input tree\n";
235  bytes = mChain->GetEntry(mEventCounter++);
236  nCycles++;
237  LOG_WARN << "Not input has been found for: " << nCycles << " times" << endm;
238  if(nCycles >= 10) {
239  LOG_ERROR << "Terminating StPicoDstReader::ReadProcess(Long64_t) after "
240  << nCycles << " times!" << endm;
241  mStatusRead = false;
242  break;
243  }
244  }
245  return mStatusRead;
246 }
Bool_t readPicoEvent(Long64_t iEvent)
Read next event in the chain.
Allows to read picoDst file(s)
void SetStatus(const Char_t *branchNameRegex, Int_t enable)
Set enable/disable branch matching when reading picoDst.
static const char * picoArrayNames[NAllPicoArrays]
Names of the TBranches in the TTree/File.
Definition: StPicoArrays.h:23
static int picoArraySizes[NAllPicoArrays]
Maximum sizes of the TClonesArrays.
Definition: StPicoArrays.h:29
static const char * picoArrayTypes[NAllPicoArrays]
Names of the classes, the TClonesArrays are arrays of this type.
Definition: StPicoArrays.h:26
static void set(TClonesArray **)
Set the pointers to the TClonesArrays.
Definition: StPicoDst.cxx:34
Main class that keeps TClonesArrays with main classes.
Definition: StPicoDst.h:40
void Init()
Calls openRead()
~StPicoDstReader()
Destructor.
void Finish()
Close files and finilize.