CMS 3D CMS Logo

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