12 #include "TIterator.h" 25 TFile *hdl = TFile::Open(fname.c_str(),
"read");
28 std::cout <<
"ERR Could not open file " << fname.c_str() << std::endl;
37 TList *keylist = hdl->GetListOfKeys();
38 TIterator *iter = keylist->MakeIterator();
40 while ((key = (TKey *)iter->Next())) {
41 TObject *
obj = hdl->Get(key->GetName());
51 TTree *
tree = (TTree *)hdl->Get(trname.c_str());
53 return tree->GetEntries();
55 std::cout <<
"ERR cannot find a TTree named \"" << trname <<
"\"" << std::endl;
61 void addBranchSizes(TBranch *
branch, Long64_t &
size) {
62 size += branch->GetTotalSize();
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);
73 if (tree !=
nullptr) {
74 Long64_t nB = tree->GetListOfBranches()->GetEntries();
75 for (Long64_t
i = 0;
i < nB; ++
i) {
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;
89 Long64_t nB = tr->GetListOfBranches()->GetEntries();
90 for (Long64_t
i = 0;
i < nB; ++
i) {
91 tr->GetListOfBranches()->At(
i)->Print();
102 uuidTree->GetEntry(0);
110 uuidTree->GetEntry(0);
117 TTree *metaDataTree) {
122 fndx->SetAddress(&findexPtr);
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";
133 std::cout <<
"\nFileFormatVersion = " << fileFormatVersion <<
". ";
134 if (fileFormatVersion.fastCopyPossible())
135 std::cout <<
"This version supports fast copy\n";
137 std::cout <<
"This version does not support fast copy\n";
139 if (fileIndex.allEventsInEntryOrder()) {
140 std::cout <<
"Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
142 std::cout <<
"Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
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";
149 std::cout <<
"Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
151 std::cout <<
"(Note that other factors can prevent fast copy from occurring)\n\n";
156 TTree *metaDataTree) {
161 fndx->SetAddress(&findexPtr);
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";
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);
177 std::cout <<
"Failed to find " << kEventAuxiliaryBranchName
178 <<
" branch in Events TTree. Something is wrong with this file." << std::endl;
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)
195 std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
206 eventAuxBranch->GetEntry(it.entry());
207 eventNum = eventAuxiliary.
id().
event();
212 std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() <<
" " << type << std::endl;
215 std::cout <<
"\nFileFormatVersion = " << fileFormatVersion <<
". ";
216 if (fileFormatVersion.fastCopyPossible())
217 std::cout <<
"This version supports fast copy\n";
219 std::cout <<
"This version does not support fast copy\n";
222 std::cout <<
"Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
224 std::cout <<
"Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
234 std::cout <<
"(Note that other factors can prevent fast copy from occurring)\n\n";
239 assert(
nullptr != metaDataTree);
245 fft->SetAddress(&fftPtr);
257 TTree *metaDataTree) {
262 fndx->SetAddress(&findexPtr);
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";
272 << std::setw(15) <<
"Run" << std::setw(15) <<
"Lumi" << std::setw(15) <<
"# Events" 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;
282 if (runID != it->run_ || lumiID != it->lumi_) {
285 std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents <<
"\n";
295 std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents <<
"\n";
303 TTree *metaDataTree) {
308 fndx->SetAddress(&findexPtr);
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";
317 << std::setw(15) <<
"Run" << std::setw(15) <<
"Lumi" << std::setw(15) <<
"# Events" 321 unsigned long runID = 0;
322 unsigned long lumiID = 0;
333 if (runID != it.run() || lumiID != it.lumi()) {
336 std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents <<
"\n";
352 std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents <<
"\n";
359 assert(
nullptr != metaDataTree);
365 fft->SetAddress(&fftPtr);
EventNumber_t event() const
void longBranchPrint(TTree *tr)
void printUuids(TTree *uuidTree)
unsigned long long EventNumber_t
std::string const & fileFormatVersionBranchName()
IndexIntoFileItr begin(SortOrder sortOrder) const
static void preIndexIntoFilePrintEventsInLumis(TFile *, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
std::string const & fileIndexBranchName()
std::string const & fid() const
std::string const & indexIntoFileBranchName()
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)
std::string getUuid(TTree *uuidTree)
const_iterator end() const
bool iterationWillBeInEntryOrder(SortOrder sortOrder) const
Used to determine whether or not to disable fast cloning.
std::string const & metaDataTreeName()
void printTrees(TFile *hdl)
void printBranchNames(TTree *tree)
void printEventLists(TFile *tfl)
void printEventsInLumis(TFile *tfl)
EventID const & id() const
std::string const & eventTreeName()
static void preIndexIntoFilePrintEventLists(TFile *, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
const_iterator begin() const
TFile * openFileHdl(std::string const &fname)
std::string const & fileIdentifierBranchName()
static void postIndexIntoFilePrintEventLists(TFile *tfl, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)
static void postIndexIntoFilePrintEventsInLumis(TFile *tfl, FileFormatVersion const &fileFormatVersion, TTree *metaDataTree)