CMS 3D CMS Logo

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