Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:08

0001 #include "IOPool/Common/bin/CollUtil.h"
0002 
0003 #include "DataFormats/Provenance/interface/BranchType.h"
0004 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
0005 #include "DataFormats/Provenance/interface/FileFormatVersion.h"
0006 #include "DataFormats/Provenance/interface/FileID.h"
0007 #include "DataFormats/Provenance/interface/FileIndex.h"
0008 #include "DataFormats/Provenance/interface/IndexIntoFile.h"
0009 
0010 #include "TBasket.h"
0011 #include "TBranch.h"
0012 #include "TFile.h"
0013 #include "TIterator.h"
0014 #include "TKey.h"
0015 #include "TList.h"
0016 #include "TObject.h"
0017 #include "TTree.h"
0018 
0019 #include <iomanip>
0020 #include <iostream>
0021 #include <ranges>
0022 
0023 namespace edm {
0024 
0025   // Get a file handler
0026   TFile *openFileHdl(std::string const &fname) {
0027     TFile *hdl = TFile::Open(fname.c_str(), "read");
0028 
0029     if (nullptr == hdl) {
0030       std::cout << "ERR Could not open file " << fname.c_str() << std::endl;
0031       exit(1);
0032     }
0033     return hdl;
0034   }
0035 
0036   // Print every tree in a file
0037   void printTrees(TFile *hdl) {
0038     hdl->ls();
0039     TList *keylist = hdl->GetListOfKeys();
0040     TIterator *iter = keylist->MakeIterator();
0041     TKey *key = nullptr;
0042     while ((key = (TKey *)iter->Next())) {
0043       TObject *obj = hdl->Get(key->GetName());
0044       if (obj->IsA() == TTree::Class()) {
0045         obj->Print();
0046       }
0047     }
0048     return;
0049   }
0050 
0051   // number of entries in a tree
0052   Long64_t numEntries(TFile *hdl, std::string const &trname) {
0053     TTree *tree = (TTree *)hdl->Get(trname.c_str());
0054     if (tree) {
0055       return tree->GetEntries();
0056     } else {
0057       std::cout << "ERR cannot find a TTree named \"" << trname << "\"" << std::endl;
0058       return -1;
0059     }
0060   }
0061 
0062   namespace {
0063     void addBranchSizes(TBranch *branch, Long64_t &size) {
0064       size += branch->GetTotalSize();  // Includes size of branch metadata
0065       // Now recurse through any subbranches.
0066       Long64_t nB = branch->GetListOfBranches()->GetEntries();
0067       for (Long64_t i = 0; i < nB; ++i) {
0068         TBranch *btemp = (TBranch *)branch->GetListOfBranches()->At(i);
0069         addBranchSizes(btemp, size);
0070       }
0071     }
0072     Long64_t branchCompressedSizes(TBranch *branch) { return branch->GetZipBytes("*"); }
0073 
0074   }  // namespace
0075 
0076   void printBranchNames(TTree *tree) {
0077     if (tree != nullptr) {
0078       Long64_t nB = tree->GetListOfBranches()->GetEntries();
0079       for (Long64_t i = 0; i < nB; ++i) {
0080         Long64_t size = 0LL;
0081         TBranch *btemp = (TBranch *)tree->GetListOfBranches()->At(i);
0082         addBranchSizes(btemp, size);
0083         std::cout << "Branch " << i << " of " << tree->GetName() << " tree: " << btemp->GetName()
0084                   << " Total size = " << size << " Compressed size = " << branchCompressedSizes(btemp) << std::endl;
0085       }
0086     } else {
0087       std::cout << "Missing Events tree?\n";
0088     }
0089   }
0090 
0091   void longBranchPrint(TTree *tr) {
0092     if (tr != nullptr) {
0093       Long64_t nB = tr->GetListOfBranches()->GetEntries();
0094       for (Long64_t i = 0; i < nB; ++i) {
0095         tr->GetListOfBranches()->At(i)->Print();
0096       }
0097     } else {
0098       std::cout << "Missing Events tree?\n";
0099     }
0100   }
0101 
0102   namespace {
0103     class BranchBasketBytes {
0104     public:
0105       BranchBasketBytes(TBranch const *branch)
0106           : basketFirstEntry_(branch->GetBasketEntry()),
0107             basketBytes_(branch->GetBasketBytes()),
0108             branchName_(branch->GetName()),
0109             maxBaskets_(branch->GetMaxBaskets()) {}
0110 
0111       bool isAlignedWithClusterBoundaries() const { return isAligned_; }
0112 
0113       std::string_view name() const { return branchName_; }
0114 
0115       // Processes "next cluster" for the branch, calculating the
0116       // number of bytes and baskets in the cluster
0117       //
0118       // @param[in] clusterBegin        Begin entry number for the cluster
0119       // @param[in] clusterEnd          End entry number (exclusive) for the cluster
0120       // @param[out] nonAlignedBranches Branch name is added to the set if the basket boundary
0121       //                                does not align with cluster boundary
0122       //
0123       // @return Tuple of the number of bytes and baskets in the cluster
0124       std::tuple<Long64_t, unsigned> bytesInNextCluster(Long64_t clusterBegin,
0125                                                         Long64_t clusterEnd,
0126                                                         std::set<std::string_view> &nonAlignedBranches) {
0127         if (basketFirstEntry_[iBasket_] != clusterBegin) {
0128           std::cout << "Branch " << branchName_ << " iBasket " << iBasket_ << " begin entry "
0129                     << basketFirstEntry_[iBasket_] << " does not align with cluster boundary, expected " << clusterBegin
0130                     << std::endl;
0131           exit(1);
0132         }
0133 
0134         Long64_t bytes = 0;
0135         unsigned nbaskets = 0;
0136         for (; iBasket_ < maxBaskets_ and basketFirstEntry_[iBasket_] < clusterEnd; ++iBasket_) {
0137           bytes += basketBytes_[iBasket_];
0138           ++nbaskets;
0139         }
0140         if (basketFirstEntry_[iBasket_] != clusterEnd) {
0141           nonAlignedBranches.insert(branchName_);
0142           isAligned_ = false;
0143           return std::tuple(0, 0);
0144         }
0145         return std::tuple(bytes, nbaskets);
0146       }
0147 
0148     private:
0149       Long64_t const *basketFirstEntry_;
0150       Int_t const *basketBytes_;
0151       std::string_view branchName_;
0152       Int_t maxBaskets_;
0153       Long64_t iBasket_ = 0;
0154       bool isAligned_ = true;
0155     };
0156 
0157     std::vector<BranchBasketBytes> makeBranchBasketBytes(TBranch *branch) {
0158       std::vector<BranchBasketBytes> ret;
0159 
0160       TObjArray *subBranches = branch->GetListOfBranches();
0161       if (subBranches and subBranches->GetEntries() > 0) {
0162         // process sub-branches if there are any
0163         auto const nbranches = subBranches->GetEntries();
0164         for (Long64_t iBranch = 0; iBranch < nbranches; ++iBranch) {
0165           auto vec = makeBranchBasketBytes(dynamic_cast<TBranch *>(subBranches->At(iBranch)));
0166           ret.insert(ret.end(), std::make_move_iterator(vec.begin()), std::make_move_iterator(vec.end()));
0167         }
0168       } else {
0169         ret.emplace_back(branch);
0170       }
0171       return ret;
0172     }
0173 
0174     template <typename T>
0175     void processClusters(TTree *tr, T printer, const std::string &limitToBranch = "") {
0176       TTree::TClusterIterator clusterIter = tr->GetClusterIterator(0);
0177       Long64_t const nentries = tr->GetEntries();
0178 
0179       // Keep the state of each branch basket index so that we don't
0180       // have to iterate through everything on every cluster
0181       std::vector<BranchBasketBytes> processors;
0182       {
0183         TObjArray *branches = tr->GetListOfBranches();
0184         Long64_t const nbranches = branches->GetEntries();
0185         for (Long64_t iBranch = 0; iBranch < nbranches; ++iBranch) {
0186           auto branch = dynamic_cast<TBranch *>(branches->At(iBranch));
0187           if (limitToBranch.empty() or
0188               std::string_view(branch->GetName()).find(limitToBranch) != std::string_view::npos) {
0189             auto vec = makeBranchBasketBytes(branch);
0190             processors.insert(
0191                 processors.end(), std::make_move_iterator(vec.begin()), std::make_move_iterator(vec.end()));
0192           }
0193         }
0194       }
0195 
0196       printer.header(tr, processors);
0197       // Record branches whose baskets do not align with cluster boundaries
0198       std::set<std::string_view> nonAlignedBranches;
0199       {
0200         Long64_t clusterBegin;
0201         while ((clusterBegin = clusterIter()) < nentries) {
0202           Long64_t clusterEnd = clusterIter.GetNextEntry();
0203           printer.beginCluster(clusterBegin, clusterEnd);
0204           for (auto &p : processors) {
0205             if (p.isAlignedWithClusterBoundaries()) {
0206               auto const [bytes, baskets] = p.bytesInNextCluster(clusterBegin, clusterEnd, nonAlignedBranches);
0207               printer.processBranch(bytes, baskets);
0208             }
0209           }
0210           printer.endCluster();
0211         }
0212       }
0213 
0214       if (not nonAlignedBranches.empty()) {
0215         std::cout << "\nThe following branches had baskets whose entry boundaries did not align with the cluster "
0216                      "boundaries. Their baskets are excluded from the cluster size calculation above starting from the "
0217                      "first basket that did not align with a cluster boundary."
0218                   << std::endl;
0219         for (auto &name : nonAlignedBranches) {
0220           std::cout << "  " << name << std::endl;
0221         }
0222       }
0223     }
0224   }  // namespace
0225 
0226   void clusterPrint(TTree *tr) {
0227     struct ClusterPrinter {
0228       void header(TTree const *tr, std::vector<BranchBasketBytes> const &branchProcessors) const {
0229         std::cout << "Printing cluster boundaries in terms of tree entries of the tree " << tr->GetName()
0230                   << ". Note that the end boundary is exclusive." << std::endl;
0231         std::cout << std::setw(15) << "Begin" << std::setw(15) << "End" << std::setw(15) << "Entries" << std::setw(15)
0232                   << "Max baskets" << std::setw(15) << "Bytes" << std::endl;
0233       }
0234 
0235       void beginCluster(Long64_t clusterBegin, Long64_t clusterEnd) {
0236         bytes_ = 0;
0237         maxbaskets_ = 0;
0238         std::cout << std::setw(15) << clusterBegin << std::setw(15) << clusterEnd << std::setw(15)
0239                   << (clusterEnd - clusterBegin);
0240       }
0241 
0242       void processBranch(Long64_t bytes, unsigned int baskets) {
0243         bytes_ += bytes;
0244         maxbaskets_ = std::max(baskets, maxbaskets_);
0245       }
0246 
0247       void endCluster() const { std::cout << std::setw(15) << maxbaskets_ << std::setw(15) << bytes_ << std::endl; }
0248 
0249       Long64_t bytes_ = 0;
0250       unsigned int maxbaskets_ = 0;
0251     };
0252     processClusters(tr, ClusterPrinter{});
0253   }
0254 
0255   void basketPrint(TTree *tr, const std::string &branchName) {
0256     struct BasketPrinter {
0257       void header(TTree const *tr, std::vector<BranchBasketBytes> const &branchProcessors) const {
0258         std::cout << "Printing cluster boundaries in terms of tree entries of the tree " << tr->GetName()
0259                   << ". Note that the end boundary is exclusive." << std::endl;
0260         std::cout << "\nBranches for which number of baskets in each cluster are printed\n";
0261         for (int i = 0; auto const &p : branchProcessors) {
0262           std::cout << "[" << i << "] " << p.name() << std::endl;
0263           ++i;
0264         }
0265         std::cout << "\n"
0266                   << std::setw(15) << "Begin" << std::setw(15) << "End" << std::setw(15) << "Entries" << std::setw(15);
0267         for (auto i : std::views::iota(0U, branchProcessors.size())) {
0268           std::cout << std::setw(5) << (std::string("[") + std::to_string(i) + "]");
0269         }
0270         std::cout << std::endl;
0271       }
0272 
0273       void beginCluster(Long64_t clusterBegin, Long64_t clusterEnd) const {
0274         std::cout << std::setw(15) << clusterBegin << std::setw(15) << clusterEnd << std::setw(15)
0275                   << (clusterEnd - clusterBegin);
0276       }
0277 
0278       void processBranch(Long64_t bytes, unsigned int baskets) const { std::cout << std::setw(5) << baskets; }
0279 
0280       void endCluster() const { std::cout << std::endl; }
0281     };
0282     processClusters(tr, BasketPrinter{}, branchName);
0283   }
0284 
0285   std::string getUuid(TTree *uuidTree) {
0286     FileID fid;
0287     FileID *fidPtr = &fid;
0288     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
0289     uuidTree->GetEntry(0);
0290     return fid.fid();
0291   }
0292 
0293   void printUuids(TTree *uuidTree) {
0294     FileID fid;
0295     FileID *fidPtr = &fid;
0296     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
0297     uuidTree->GetEntry(0);
0298 
0299     std::cout << "UUID: " << fid.fid() << std::endl;
0300   }
0301 
0302   static void preIndexIntoFilePrintEventLists(TFile *,
0303                                               FileFormatVersion const &fileFormatVersion,
0304                                               TTree *metaDataTree) {
0305     FileIndex fileIndex;
0306     FileIndex *findexPtr = &fileIndex;
0307     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != nullptr) {
0308       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
0309       fndx->SetAddress(&findexPtr);
0310       fndx->GetEntry(0);
0311     } else {
0312       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
0313                    "this indicates a problem with the file.  This condition should be expected with\n"
0314                    "files created with earlier releases and printout of the event list will fail.\n";
0315       return;
0316     }
0317 
0318     std::cout << "\n" << fileIndex;
0319 
0320     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
0321     if (fileFormatVersion.fastCopyPossible())
0322       std::cout << "This version supports fast copy\n";
0323     else
0324       std::cout << "This version does not support fast copy\n";
0325 
0326     if (fileIndex.allEventsInEntryOrder()) {
0327       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
0328     } else {
0329       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
0330     }
0331 
0332     fileIndex.sortBy_Run_Lumi_EventEntry();
0333     if (fileIndex.allEventsInEntryOrder()) {
0334       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
0335     } else {
0336       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
0337     }
0338     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
0339   }
0340 
0341   static void postIndexIntoFilePrintEventLists(TFile *tfl,
0342                                                FileFormatVersion const &fileFormatVersion,
0343                                                TTree *metaDataTree) {
0344     IndexIntoFile indexIntoFile;
0345     IndexIntoFile *findexPtr = &indexIntoFile;
0346     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != nullptr) {
0347       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
0348       fndx->SetAddress(&findexPtr);
0349       fndx->GetEntry(0);
0350     } else {
0351       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
0352                    "this indicates a problem with the file.  This condition should be expected with\n"
0353                    "files created with earlier releases and printout of the event list will fail.\n";
0354       return;
0355     }
0356     //need to read event # from the EventAuxiliary branch
0357     TTree *eventsTree = dynamic_cast<TTree *>(tfl->Get(poolNames::eventTreeName().c_str()));
0358     TBranch *eventAuxBranch = nullptr;
0359     assert(nullptr != eventsTree);
0360     char const *const kEventAuxiliaryBranchName = "EventAuxiliary";
0361     if (eventsTree->FindBranch(kEventAuxiliaryBranchName) != nullptr) {
0362       eventAuxBranch = eventsTree->GetBranch(kEventAuxiliaryBranchName);
0363     } else {
0364       std::cout << "Failed to find " << kEventAuxiliaryBranchName
0365                 << " branch in Events TTree.  Something is wrong with this file." << std::endl;
0366       return;
0367     }
0368     EventAuxiliary eventAuxiliary;
0369     EventAuxiliary *eAPtr = &eventAuxiliary;
0370     eventAuxBranch->SetAddress(&eAPtr);
0371     std::cout << "\nPrinting IndexIntoFile contents.  This includes a list of all Runs, LuminosityBlocks\n"
0372               << "and Events stored in the root file.\n\n";
0373     std::cout << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "Event" << std::setw(15)
0374               << "TTree Entry"
0375               << "\n";
0376 
0377     for (IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
0378                                          itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
0379          it != itEnd;
0380          ++it) {
0381       IndexIntoFile::EntryType t = it.getEntryType();
0382       std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
0383       EventNumber_t eventNum = 0;
0384       std::string type;
0385       switch (t) {
0386         case IndexIntoFile::kRun:
0387           type = "(Run)";
0388           break;
0389         case IndexIntoFile::kLumi:
0390           type = "(Lumi)";
0391           break;
0392         case IndexIntoFile::kEvent:
0393           eventAuxBranch->GetEntry(it.entry());
0394           eventNum = eventAuxiliary.id().event();
0395           break;
0396         default:
0397           break;
0398       }
0399       std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() << " " << type << std::endl;
0400     }
0401 
0402     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
0403     if (fileFormatVersion.fastCopyPossible())
0404       std::cout << "This version supports fast copy\n";
0405     else
0406       std::cout << "This version does not support fast copy\n";
0407 
0408     if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::firstAppearanceOrder)) {
0409       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
0410     } else {
0411       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
0412     }
0413 
0414     // This will not work unless the other nonpersistent parts of the Index are filled first
0415     // I did not have time to implement this yet.
0416     // if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::numericalOrder)) {
0417     //   std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
0418     // } else {
0419     //   std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
0420     // }
0421     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
0422   }
0423 
0424   void printEventLists(TFile *tfl) {
0425     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
0426     assert(nullptr != metaDataTree);
0427 
0428     FileFormatVersion fileFormatVersion;
0429     FileFormatVersion *fftPtr = &fileFormatVersion;
0430     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != nullptr) {
0431       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
0432       fft->SetAddress(&fftPtr);
0433       fft->GetEntry(0);
0434     }
0435     if (fileFormatVersion.hasIndexIntoFile()) {
0436       postIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
0437     } else {
0438       preIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
0439     }
0440   }
0441 
0442   static void preIndexIntoFilePrintEventsInLumis(TFile *,
0443                                                  FileFormatVersion const &fileFormatVersion,
0444                                                  TTree *metaDataTree) {
0445     FileIndex fileIndex;
0446     FileIndex *findexPtr = &fileIndex;
0447     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != nullptr) {
0448       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
0449       fndx->SetAddress(&findexPtr);
0450       fndx->GetEntry(0);
0451     } else {
0452       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
0453                    "this indicates a problem with the file.  This condition should be expected with\n"
0454                    "files created with earlier releases and printout of the event list will fail.\n";
0455       return;
0456     }
0457 
0458     std::cout << "\n"
0459               << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "# Events"
0460               << "\n";
0461     unsigned long nEvents = 0;
0462     unsigned long runID = 0;
0463     unsigned long lumiID = 0;
0464     for (std::vector<FileIndex::Element>::const_iterator it = fileIndex.begin(), itEnd = fileIndex.end(); it != itEnd;
0465          ++it) {
0466       if (it->getEntryType() == FileIndex::kEvent) {
0467         ++nEvents;
0468       } else if (it->getEntryType() == FileIndex::kLumi) {
0469         if (runID != it->run_ || lumiID != it->lumi_) {
0470           //print the previous one
0471           if (lumiID != 0) {
0472             std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0473           }
0474           nEvents = 0;
0475           runID = it->run_;
0476           lumiID = it->lumi_;
0477         }
0478       }
0479     }
0480     //print the last one
0481     if (lumiID != 0) {
0482       std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0483     }
0484 
0485     std::cout << "\n";
0486   }
0487 
0488   static void postIndexIntoFilePrintEventsInLumis(TFile *tfl,
0489                                                   FileFormatVersion const &fileFormatVersion,
0490                                                   TTree *metaDataTree) {
0491     IndexIntoFile indexIntoFile;
0492     IndexIntoFile *findexPtr = &indexIntoFile;
0493     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != nullptr) {
0494       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
0495       fndx->SetAddress(&findexPtr);
0496       fndx->GetEntry(0);
0497     } else {
0498       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
0499                    "this indicates a problem with the file.  This condition should be expected with\n"
0500                    "files created with earlier releases and printout of the event list will fail.\n";
0501       return;
0502     }
0503     std::cout << "\n"
0504               << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "# Events"
0505               << "\n";
0506 
0507     unsigned long nEvents = 0;
0508     unsigned long runID = 0;
0509     unsigned long lumiID = 0;
0510 
0511     for (IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
0512                                          itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
0513          it != itEnd;
0514          ++it) {
0515       IndexIntoFile::EntryType t = it.getEntryType();
0516       switch (t) {
0517         case IndexIntoFile::kRun:
0518           break;
0519         case IndexIntoFile::kLumi:
0520           if (runID != it.run() || lumiID != it.lumi()) {
0521             //print the previous one
0522             if (lumiID != 0) {
0523               std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0524             }
0525             nEvents = 0;
0526             runID = it.run();
0527             lumiID = it.lumi();
0528           }
0529           break;
0530         case IndexIntoFile::kEvent:
0531           ++nEvents;
0532           break;
0533         default:
0534           break;
0535       }
0536     }
0537     //print the last one
0538     if (lumiID != 0) {
0539       std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0540     }
0541     std::cout << "\n";
0542   }
0543 
0544   void printEventsInLumis(TFile *tfl) {
0545     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
0546     assert(nullptr != metaDataTree);
0547 
0548     FileFormatVersion fileFormatVersion;
0549     FileFormatVersion *fftPtr = &fileFormatVersion;
0550     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != nullptr) {
0551       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
0552       fft->SetAddress(&fftPtr);
0553       fft->GetEntry(0);
0554     }
0555     if (fileFormatVersion.hasIndexIntoFile()) {
0556       postIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
0557     } else {
0558       preIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
0559     }
0560   }
0561 }  // namespace edm