CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CollUtil.cc
Go to the documentation of this file.
2 
9 
10 #include "TBranch.h"
11 #include "TFile.h"
12 #include "TIterator.h"
13 #include "TKey.h"
14 #include "TList.h"
15 #include "TObject.h"
16 #include "TTree.h"
17 
18 #include <iomanip>
19 #include <iostream>
20 
21 namespace edm {
22 
23  // Get a file handler
24  TFile* openFileHdl(std::string const& fname) {
25  TFile *hdl = TFile::Open(fname.c_str(), "read");
26 
27  if (0 == hdl) {
28  std::cout << "ERR Could not open file " << fname.c_str() << std::endl;
29  exit(1);
30  }
31  return hdl;
32  }
33 
34  // Print every tree in a file
35  void printTrees(TFile *hdl) {
36  hdl->ls();
37  TList *keylist = hdl->GetListOfKeys();
38  TIterator *iter = keylist->MakeIterator();
39  TKey *key = 0;
40  while ((key = (TKey*)iter->Next())) {
41  TObject *obj = hdl->Get(key->GetName());
42  if (obj->IsA() == TTree::Class()) {
43  obj->Print();
44  }
45  }
46  return;
47  }
48 
49  // number of entries in a tree
50  Long64_t numEntries(TFile *hdl, std::string const& trname) {
51  TTree *tree = (TTree*)hdl->Get(trname.c_str());
52  if (tree) {
53  return tree->GetEntries();
54  } else {
55  std::cout << "ERR cannot find a TTree named \"" << trname << "\""
56  << std::endl;
57  return -1;
58  }
59  }
60 
61  namespace {
62  void addBranchSizes(TBranch *branch, Long64_t& size) {
63  size += branch->GetTotalSize(); // Includes size of branch metadata
64  // Now recurse through any subbranches.
65  Long64_t nB = branch->GetListOfBranches()->GetEntries();
66  for (Long64_t i = 0; i < nB; ++i) {
67  TBranch *btemp = (TBranch *)branch->GetListOfBranches()->At(i);
68  addBranchSizes(btemp, size);
69  }
70  }
71  }
72 
73  void printBranchNames(TTree *tree) {
74  if (tree != 0) {
75  Long64_t nB = tree->GetListOfBranches()->GetEntries();
76  for (Long64_t i = 0; i < nB; ++i) {
77  Long64_t size = 0LL;
78  TBranch *btemp = (TBranch *)tree->GetListOfBranches()->At(i);
79  addBranchSizes(btemp, size);
80  std::cout << "Branch " << i << " of " << tree->GetName() << " tree: " << btemp->GetName() << " Total size = " << size << std::endl;
81  }
82  } else {
83  std::cout << "Missing Events tree?\n";
84  }
85  }
86 
87  void longBranchPrint(TTree *tr) {
88  if (tr != 0) {
89  Long64_t nB = tr->GetListOfBranches()->GetEntries();
90  for (Long64_t i = 0; i < nB; ++i) {
91  tr->GetListOfBranches()->At(i)->Print();
92  }
93  } else {
94  std::cout << "Missing Events tree?\n";
95  }
96  }
97 
98  std::string getUuid(TTree *uuidTree) {
99  FileID fid;
100  FileID *fidPtr = &fid;
101  uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
102  uuidTree->GetEntry(0);
103  return fid.fid();
104  }
105 
106  void printUuids(TTree *uuidTree) {
107  FileID fid;
108  FileID *fidPtr = &fid;
109  uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
110  uuidTree->GetEntry(0);
111 
112  std::cout << "UUID: " << fid.fid() << std::endl;
113  }
114 
115  static void preIndexIntoFilePrintEventLists(TFile*, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
116  FileIndex fileIndex;
117  FileIndex *findexPtr = &fileIndex;
118  if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != 0) {
119  TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
120  fndx->SetAddress(&findexPtr);
121  fndx->GetEntry(0);
122  } else {
123  std::cout << "FileIndex not found. If this input file was created with release 1_8_0 or later\n"
124  "this indicates a problem with the file. This condition should be expected with\n"
125  "files created with earlier releases and printout of the event list will fail.\n";
126  return;
127  }
128 
129  std::cout << "\n" << fileIndex;
130 
131  std::cout << "\nFileFormatVersion = " << fileFormatVersion << ". ";
132  if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
133  else std::cout << "This version does not support fast copy\n";
134 
135  if (fileIndex.allEventsInEntryOrder()) {
136  std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
137  } else {
138  std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
139  }
140 
141  fileIndex.sortBy_Run_Lumi_EventEntry();
142  if (fileIndex.allEventsInEntryOrder()) {
143  std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
144  } else {
145  std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
146  }
147  std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
148  }
149 
150  static void postIndexIntoFilePrintEventLists(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
151  IndexIntoFile indexIntoFile;
152  IndexIntoFile *findexPtr = &indexIntoFile;
153  if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != 0) {
154  TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
155  fndx->SetAddress(&findexPtr);
156  fndx->GetEntry(0);
157  } else {
158  std::cout << "IndexIntoFile not found. If this input file was created with release 1_8_0 or later\n"
159  "this indicates a problem with the file. This condition should be expected with\n"
160  "files created with earlier releases and printout of the event list will fail.\n";
161  return;
162  }
163  //need to read event # from the EventAuxiliary branch
164  TTree* eventsTree = dynamic_cast<TTree*>(tfl->Get(poolNames::eventTreeName().c_str()));
165  TBranch* eventAuxBranch = 0;
166  assert(0 != eventsTree);
167  char const* const kEventAuxiliaryBranchName = "EventAuxiliary";
168  if(eventsTree->FindBranch(kEventAuxiliaryBranchName) != 0){
169  eventAuxBranch = eventsTree->GetBranch(kEventAuxiliaryBranchName);
170  } else {
171  std::cout << "Failed to find " << kEventAuxiliaryBranchName << " branch in Events TTree. Something is wrong with this file." << std::endl;
172  return;
173  }
174  EventAuxiliary eventAuxiliary;
175  EventAuxiliary* eAPtr = &eventAuxiliary;
176  eventAuxBranch->SetAddress(&eAPtr);
177  std::cout << "\nPrinting IndexIntoFile contents. This includes a list of all Runs, LuminosityBlocks\n"
178  << "and Events stored in the root file.\n\n";
179  std::cout << std::setw(15) << "Run"
180  << std::setw(15) << "Lumi"
181  << std::setw(15) << "Event"
182  << std::setw(15) << "TTree Entry"
183  << "\n";
184 
186  itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
187  it != itEnd; ++it) {
188  IndexIntoFile::EntryType t = it.getEntryType();
189  std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
190  EventNumber_t eventNum = 0;
192  switch(t) {
193  case IndexIntoFile::kRun:
194  type = "(Run)";
195  break;
197  type = "(Lumi)";
198  break;
200  eventAuxBranch->GetEntry(it.entry());
201  eventNum = eventAuxiliary.id().event();
202  break;
203  default:
204  break;
205  }
206  std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() << " " << type << std::endl;
207  }
208 
209  std::cout << "\nFileFormatVersion = " << fileFormatVersion << ". ";
210  if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
211  else std::cout << "This version does not support fast copy\n";
212 
214  std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
215  } else {
216  std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
217  }
218 
219  // This will not work unless the other nonpersistent parts of the Index are filled first
220  // I did not have time to implement this yet.
221  // if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::numericalOrder)) {
222  // std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
223  // } else {
224  // std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
225  // }
226  std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
227  }
228 
229  void printEventLists(TFile *tfl) {
230  TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
231  assert(0 != metaDataTree);
232 
233  FileFormatVersion fileFormatVersion;
234  FileFormatVersion *fftPtr = &fileFormatVersion;
235  if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != 0) {
236  TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
237  fft->SetAddress(&fftPtr);
238  fft->GetEntry(0);
239  }
240  if(fileFormatVersion.hasIndexIntoFile()) {
241  postIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
242  } else {
243  preIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
244  }
245  }
246 
247  static void preIndexIntoFilePrintEventsInLumis(TFile*, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
248  FileIndex fileIndex;
249  FileIndex *findexPtr = &fileIndex;
250  if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != 0) {
251  TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
252  fndx->SetAddress(&findexPtr);
253  fndx->GetEntry(0);
254  } else {
255  std::cout << "FileIndex not found. If this input file was created with release 1_8_0 or later\n"
256  "this indicates a problem with the file. This condition should be expected with\n"
257  "files created with earlier releases and printout of the event list will fail.\n";
258  return;
259  }
260 
261  std::cout <<"\n"<< std::setw(15) << "Run"
262  << std::setw(15) << "Lumi"
263  << std::setw(15) << "# Events"
264  << "\n";
265  unsigned long nEvents = 0;
266  unsigned long runID = 0;
267  unsigned long lumiID = 0;
268  for(std::vector<FileIndex::Element>::const_iterator it = fileIndex.begin(), itEnd = fileIndex.end(); it != itEnd; ++it) {
269  if(it->getEntryType() == FileIndex::kEvent) {
270  ++nEvents;
271  }
272  else if(it->getEntryType() == FileIndex::kLumi) {
273  if(runID !=it->run_ || lumiID !=it->lumi_) {
274  //print the previous one
275  if(lumiID !=0) {
276  std::cout << std::setw(15) << runID
277  << std::setw(15) << lumiID
278  << std::setw(15) << nEvents<<"\n";
279  }
280  nEvents=0;
281  runID = it->run_;
282  lumiID = it->lumi_;
283  }
284  }
285  }
286  //print the last one
287  if(lumiID !=0) {
288  std::cout << std::setw(15) << runID
289  << std::setw(15) << lumiID
290  << std::setw(15) << nEvents<<"\n";
291  }
292 
293  std::cout << "\n";
294  }
295 
296  static void postIndexIntoFilePrintEventsInLumis(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
297  IndexIntoFile indexIntoFile;
298  IndexIntoFile *findexPtr = &indexIntoFile;
299  if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != 0) {
300  TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
301  fndx->SetAddress(&findexPtr);
302  fndx->GetEntry(0);
303  } else {
304  std::cout << "IndexIntoFile not found. If this input file was created with release 1_8_0 or later\n"
305  "this indicates a problem with the file. This condition should be expected with\n"
306  "files created with earlier releases and printout of the event list will fail.\n";
307  return;
308  }
309  std::cout <<"\n"<< std::setw(15) << "Run"
310  << std::setw(15) << "Lumi"
311  << std::setw(15) << "# Events"
312  << "\n";
313 
314  unsigned long nEvents = 0;
315  unsigned long runID = 0;
316  unsigned long lumiID = 0;
317 
319  itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
320  it != itEnd; ++it) {
321  IndexIntoFile::EntryType t = it.getEntryType();
322  switch(t) {
323  case IndexIntoFile::kRun:
324  break;
326  if(runID != it.run() || lumiID != it.lumi()) {
327  //print the previous one
328  if(lumiID !=0) {
329  std::cout << std::setw(15) << runID
330  << std::setw(15) << lumiID
331  << std::setw(15) << nEvents<<"\n";
332  }
333  nEvents=0;
334  runID = it.run();
335  lumiID = it.lumi();
336  }
337  break;
339  ++nEvents;
340  break;
341  default:
342  break;
343  }
344  }
345  //print the last one
346  if(lumiID !=0) {
347  std::cout << std::setw(15) << runID
348  << std::setw(15) << lumiID
349  << std::setw(15) << nEvents<<"\n";
350  }
351  std::cout << "\n";
352  }
353 
354  void printEventsInLumis(TFile *tfl) {
355  TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
356  assert(0 != metaDataTree);
357 
358  FileFormatVersion fileFormatVersion;
359  FileFormatVersion *fftPtr = &fileFormatVersion;
360  if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != 0) {
361  TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
362  fft->SetAddress(&fftPtr);
363  fft->GetEntry(0);
364  }
365  if(fileFormatVersion.hasIndexIntoFile()) {
366  postIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
367  } else {
368  preIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
369  }
370  }
371 }
type
Definition: HCALResponse.h:21
EventNumber_t event() const
Definition: EventID.h:41
void longBranchPrint(TTree *tr)
Definition: CollUtil.cc:87
int i
Definition: DBlmapReader.cc:9
void printUuids(TTree *uuidTree)
Definition: CollUtil.cc:106
bool hasIndexIntoFile() const
assert(m_qm.get())
unsigned long long EventNumber_t
std::string const & fileFormatVersionBranchName()
Definition: BranchType.cc:218
IndexIntoFileItr begin(SortOrder sortOrder) const
TFile * openFileHdl(std::string const &fname)
Definition: CollUtil.cc:24
static void preIndexIntoFilePrintEventsInLumis(TFile *, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
Definition: CollUtil.cc:247
std::string const & fileIndexBranchName()
Definition: BranchType.cc:228
std::string const & fid() const
Definition: FileID.h:19
std::string const & indexIntoFileBranchName()
Definition: BranchType.cc:233
IndexIntoFileItr end(SortOrder sortOrder) const
Used to end an iteration over the Runs, Lumis, and Events in a file.
Long64_t numEntries(TFile *hdl, std::string const &trname)
Definition: CollUtil.cc:50
std::string getUuid(TTree *uuidTree)
Definition: CollUtil.cc:98
const_iterator end() const
Definition: FileIndex.h:100
bool iterationWillBeInEntryOrder(SortOrder sortOrder) const
Used to determine whether or not to disable fast cloning.
std::string const & metaDataTreeName()
Definition: BranchType.cc:168
void printTrees(TFile *hdl)
Definition: CollUtil.cc:35
void printBranchNames(TTree *tree)
Definition: CollUtil.cc:73
void printEventLists(TFile *tfl)
Definition: CollUtil.cc:229
string fname
main script
void printEventsInLumis(TFile *tfl)
Definition: CollUtil.cc:354
EventID const & id() const
std::string const & eventTreeName()
Definition: BranchType.cc:260
static void preIndexIntoFilePrintEventLists(TFile *, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
Definition: CollUtil.cc:115
tuple cout
Definition: gather_cfg.py:121
const_iterator begin() const
Definition: FileIndex.h:98
UInt_t nEvents
Definition: hcalCalib.cc:42
list fid
Definition: NewTree.py:51
std::string const & fileIdentifierBranchName()
Definition: BranchType.cc:223
tuple size
Write out results.
static void postIndexIntoFilePrintEventLists(TFile *tfl, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
Definition: CollUtil.cc:150
static void postIndexIntoFilePrintEventsInLumis(TFile *tfl, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
Definition: CollUtil.cc:296