Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:21:11

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 "TBranch.h"
0011 #include "TFile.h"
0012 #include "TIterator.h"
0013 #include "TKey.h"
0014 #include "TList.h"
0015 #include "TObject.h"
0016 #include "TTree.h"
0017 
0018 #include <iomanip>
0019 #include <iostream>
0020 
0021 namespace edm {
0022 
0023   // Get a file handler
0024   TFile *openFileHdl(std::string const &fname) {
0025     TFile *hdl = TFile::Open(fname.c_str(), "read");
0026 
0027     if (nullptr == hdl) {
0028       std::cout << "ERR Could not open file " << fname.c_str() << std::endl;
0029       exit(1);
0030     }
0031     return hdl;
0032   }
0033 
0034   // Print every tree in a file
0035   void printTrees(TFile *hdl) {
0036     hdl->ls();
0037     TList *keylist = hdl->GetListOfKeys();
0038     TIterator *iter = keylist->MakeIterator();
0039     TKey *key = nullptr;
0040     while ((key = (TKey *)iter->Next())) {
0041       TObject *obj = hdl->Get(key->GetName());
0042       if (obj->IsA() == TTree::Class()) {
0043         obj->Print();
0044       }
0045     }
0046     return;
0047   }
0048 
0049   // number of entries in a tree
0050   Long64_t numEntries(TFile *hdl, std::string const &trname) {
0051     TTree *tree = (TTree *)hdl->Get(trname.c_str());
0052     if (tree) {
0053       return tree->GetEntries();
0054     } else {
0055       std::cout << "ERR cannot find a TTree named \"" << trname << "\"" << std::endl;
0056       return -1;
0057     }
0058   }
0059 
0060   namespace {
0061     void addBranchSizes(TBranch *branch, Long64_t &size) {
0062       size += branch->GetTotalSize();  // Includes size of branch metadata
0063       // Now recurse through any subbranches.
0064       Long64_t nB = branch->GetListOfBranches()->GetEntries();
0065       for (Long64_t i = 0; i < nB; ++i) {
0066         TBranch *btemp = (TBranch *)branch->GetListOfBranches()->At(i);
0067         addBranchSizes(btemp, size);
0068       }
0069     }
0070   }  // namespace
0071 
0072   void printBranchNames(TTree *tree) {
0073     if (tree != nullptr) {
0074       Long64_t nB = tree->GetListOfBranches()->GetEntries();
0075       for (Long64_t i = 0; i < nB; ++i) {
0076         Long64_t size = 0LL;
0077         TBranch *btemp = (TBranch *)tree->GetListOfBranches()->At(i);
0078         addBranchSizes(btemp, size);
0079         std::cout << "Branch " << i << " of " << tree->GetName() << " tree: " << btemp->GetName()
0080                   << " Total size = " << size << std::endl;
0081       }
0082     } else {
0083       std::cout << "Missing Events tree?\n";
0084     }
0085   }
0086 
0087   void longBranchPrint(TTree *tr) {
0088     if (tr != nullptr) {
0089       Long64_t nB = tr->GetListOfBranches()->GetEntries();
0090       for (Long64_t i = 0; i < nB; ++i) {
0091         tr->GetListOfBranches()->At(i)->Print();
0092       }
0093     } else {
0094       std::cout << "Missing Events tree?\n";
0095     }
0096   }
0097 
0098   std::string getUuid(TTree *uuidTree) {
0099     FileID fid;
0100     FileID *fidPtr = &fid;
0101     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
0102     uuidTree->GetEntry(0);
0103     return fid.fid();
0104   }
0105 
0106   void printUuids(TTree *uuidTree) {
0107     FileID fid;
0108     FileID *fidPtr = &fid;
0109     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
0110     uuidTree->GetEntry(0);
0111 
0112     std::cout << "UUID: " << fid.fid() << std::endl;
0113   }
0114 
0115   static void preIndexIntoFilePrintEventLists(TFile *,
0116                                               FileFormatVersion const &fileFormatVersion,
0117                                               TTree *metaDataTree) {
0118     FileIndex fileIndex;
0119     FileIndex *findexPtr = &fileIndex;
0120     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != nullptr) {
0121       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
0122       fndx->SetAddress(&findexPtr);
0123       fndx->GetEntry(0);
0124     } else {
0125       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
0126                    "this indicates a problem with the file.  This condition should be expected with\n"
0127                    "files created with earlier releases and printout of the event list will fail.\n";
0128       return;
0129     }
0130 
0131     std::cout << "\n" << fileIndex;
0132 
0133     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
0134     if (fileFormatVersion.fastCopyPossible())
0135       std::cout << "This version supports fast copy\n";
0136     else
0137       std::cout << "This version does not support fast copy\n";
0138 
0139     if (fileIndex.allEventsInEntryOrder()) {
0140       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
0141     } else {
0142       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
0143     }
0144 
0145     fileIndex.sortBy_Run_Lumi_EventEntry();
0146     if (fileIndex.allEventsInEntryOrder()) {
0147       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
0148     } else {
0149       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
0150     }
0151     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
0152   }
0153 
0154   static void postIndexIntoFilePrintEventLists(TFile *tfl,
0155                                                FileFormatVersion const &fileFormatVersion,
0156                                                TTree *metaDataTree) {
0157     IndexIntoFile indexIntoFile;
0158     IndexIntoFile *findexPtr = &indexIntoFile;
0159     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != nullptr) {
0160       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
0161       fndx->SetAddress(&findexPtr);
0162       fndx->GetEntry(0);
0163     } else {
0164       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
0165                    "this indicates a problem with the file.  This condition should be expected with\n"
0166                    "files created with earlier releases and printout of the event list will fail.\n";
0167       return;
0168     }
0169     //need to read event # from the EventAuxiliary branch
0170     TTree *eventsTree = dynamic_cast<TTree *>(tfl->Get(poolNames::eventTreeName().c_str()));
0171     TBranch *eventAuxBranch = nullptr;
0172     assert(nullptr != eventsTree);
0173     char const *const kEventAuxiliaryBranchName = "EventAuxiliary";
0174     if (eventsTree->FindBranch(kEventAuxiliaryBranchName) != nullptr) {
0175       eventAuxBranch = eventsTree->GetBranch(kEventAuxiliaryBranchName);
0176     } else {
0177       std::cout << "Failed to find " << kEventAuxiliaryBranchName
0178                 << " branch in Events TTree.  Something is wrong with this file." << std::endl;
0179       return;
0180     }
0181     EventAuxiliary eventAuxiliary;
0182     EventAuxiliary *eAPtr = &eventAuxiliary;
0183     eventAuxBranch->SetAddress(&eAPtr);
0184     std::cout << "\nPrinting IndexIntoFile contents.  This includes a list of all Runs, LuminosityBlocks\n"
0185               << "and Events stored in the root file.\n\n";
0186     std::cout << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "Event" << std::setw(15)
0187               << "TTree Entry"
0188               << "\n";
0189 
0190     for (IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
0191                                          itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
0192          it != itEnd;
0193          ++it) {
0194       IndexIntoFile::EntryType t = it.getEntryType();
0195       std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
0196       EventNumber_t eventNum = 0;
0197       std::string type;
0198       switch (t) {
0199         case IndexIntoFile::kRun:
0200           type = "(Run)";
0201           break;
0202         case IndexIntoFile::kLumi:
0203           type = "(Lumi)";
0204           break;
0205         case IndexIntoFile::kEvent:
0206           eventAuxBranch->GetEntry(it.entry());
0207           eventNum = eventAuxiliary.id().event();
0208           break;
0209         default:
0210           break;
0211       }
0212       std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() << " " << type << std::endl;
0213     }
0214 
0215     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
0216     if (fileFormatVersion.fastCopyPossible())
0217       std::cout << "This version supports fast copy\n";
0218     else
0219       std::cout << "This version does not support fast copy\n";
0220 
0221     if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::firstAppearanceOrder)) {
0222       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
0223     } else {
0224       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
0225     }
0226 
0227     // This will not work unless the other nonpersistent parts of the Index are filled first
0228     // I did not have time to implement this yet.
0229     // if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::numericalOrder)) {
0230     //   std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
0231     // } else {
0232     //   std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
0233     // }
0234     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
0235   }
0236 
0237   void printEventLists(TFile *tfl) {
0238     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
0239     assert(nullptr != metaDataTree);
0240 
0241     FileFormatVersion fileFormatVersion;
0242     FileFormatVersion *fftPtr = &fileFormatVersion;
0243     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != nullptr) {
0244       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
0245       fft->SetAddress(&fftPtr);
0246       fft->GetEntry(0);
0247     }
0248     if (fileFormatVersion.hasIndexIntoFile()) {
0249       postIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
0250     } else {
0251       preIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
0252     }
0253   }
0254 
0255   static void preIndexIntoFilePrintEventsInLumis(TFile *,
0256                                                  FileFormatVersion const &fileFormatVersion,
0257                                                  TTree *metaDataTree) {
0258     FileIndex fileIndex;
0259     FileIndex *findexPtr = &fileIndex;
0260     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != nullptr) {
0261       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
0262       fndx->SetAddress(&findexPtr);
0263       fndx->GetEntry(0);
0264     } else {
0265       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
0266                    "this indicates a problem with the file.  This condition should be expected with\n"
0267                    "files created with earlier releases and printout of the event list will fail.\n";
0268       return;
0269     }
0270 
0271     std::cout << "\n"
0272               << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "# Events"
0273               << "\n";
0274     unsigned long nEvents = 0;
0275     unsigned long runID = 0;
0276     unsigned long lumiID = 0;
0277     for (std::vector<FileIndex::Element>::const_iterator it = fileIndex.begin(), itEnd = fileIndex.end(); it != itEnd;
0278          ++it) {
0279       if (it->getEntryType() == FileIndex::kEvent) {
0280         ++nEvents;
0281       } else if (it->getEntryType() == FileIndex::kLumi) {
0282         if (runID != it->run_ || lumiID != it->lumi_) {
0283           //print the previous one
0284           if (lumiID != 0) {
0285             std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0286           }
0287           nEvents = 0;
0288           runID = it->run_;
0289           lumiID = it->lumi_;
0290         }
0291       }
0292     }
0293     //print the last one
0294     if (lumiID != 0) {
0295       std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0296     }
0297 
0298     std::cout << "\n";
0299   }
0300 
0301   static void postIndexIntoFilePrintEventsInLumis(TFile *tfl,
0302                                                   FileFormatVersion const &fileFormatVersion,
0303                                                   TTree *metaDataTree) {
0304     IndexIntoFile indexIntoFile;
0305     IndexIntoFile *findexPtr = &indexIntoFile;
0306     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != nullptr) {
0307       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
0308       fndx->SetAddress(&findexPtr);
0309       fndx->GetEntry(0);
0310     } else {
0311       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
0312                    "this indicates a problem with the file.  This condition should be expected with\n"
0313                    "files created with earlier releases and printout of the event list will fail.\n";
0314       return;
0315     }
0316     std::cout << "\n"
0317               << std::setw(15) << "Run" << std::setw(15) << "Lumi" << std::setw(15) << "# Events"
0318               << "\n";
0319 
0320     unsigned long nEvents = 0;
0321     unsigned long runID = 0;
0322     unsigned long lumiID = 0;
0323 
0324     for (IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
0325                                          itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
0326          it != itEnd;
0327          ++it) {
0328       IndexIntoFile::EntryType t = it.getEntryType();
0329       switch (t) {
0330         case IndexIntoFile::kRun:
0331           break;
0332         case IndexIntoFile::kLumi:
0333           if (runID != it.run() || lumiID != it.lumi()) {
0334             //print the previous one
0335             if (lumiID != 0) {
0336               std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0337             }
0338             nEvents = 0;
0339             runID = it.run();
0340             lumiID = it.lumi();
0341           }
0342           break;
0343         case IndexIntoFile::kEvent:
0344           ++nEvents;
0345           break;
0346         default:
0347           break;
0348       }
0349     }
0350     //print the last one
0351     if (lumiID != 0) {
0352       std::cout << std::setw(15) << runID << std::setw(15) << lumiID << std::setw(15) << nEvents << "\n";
0353     }
0354     std::cout << "\n";
0355   }
0356 
0357   void printEventsInLumis(TFile *tfl) {
0358     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
0359     assert(nullptr != metaDataTree);
0360 
0361     FileFormatVersion fileFormatVersion;
0362     FileFormatVersion *fftPtr = &fileFormatVersion;
0363     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != nullptr) {
0364       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
0365       fft->SetAddress(&fftPtr);
0366       fft->GetEntry(0);
0367     }
0368     if (fileFormatVersion.hasIndexIntoFile()) {
0369       postIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
0370     } else {
0371       preIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
0372     }
0373   }
0374 }  // namespace edm