Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:04

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