CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/IOPool/Common/bin/CollUtil.cc

Go to the documentation of this file.
00001 #include "IOPool/Common/bin/CollUtil.h"
00002 
00003 #include "DataFormats/Provenance/interface/BranchType.h"
00004 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
00005 #include "DataFormats/Provenance/interface/FileFormatVersion.h"
00006 #include "DataFormats/Provenance/interface/FileID.h"
00007 #include "DataFormats/Provenance/interface/FileIndex.h"
00008 #include "DataFormats/Provenance/interface/IndexIntoFile.h"
00009 
00010 #include "TBranch.h"
00011 #include "TFile.h"
00012 #include "TIterator.h"
00013 #include "TKey.h"
00014 #include "TList.h"
00015 #include "TObject.h"
00016 #include "TTree.h"
00017 
00018 #include <iomanip>
00019 #include <iostream>
00020 
00021 namespace edm {
00022 
00023   // Get a file handler
00024   TFile* openFileHdl(std::string const& fname) {
00025     TFile *hdl = TFile::Open(fname.c_str(), "read");
00026 
00027     if (0 == hdl) {
00028       std::cout << "ERR Could not open file " << fname.c_str() << std::endl;
00029       exit(1);
00030     }
00031     return hdl;
00032   }
00033 
00034   // Print every tree in a file
00035   void printTrees(TFile *hdl) {
00036     hdl->ls();
00037     TList *keylist = hdl->GetListOfKeys();
00038     TIterator *iter = keylist->MakeIterator();
00039     TKey *key = 0;
00040     while ((key = (TKey*)iter->Next())) {
00041       TObject *obj = hdl->Get(key->GetName());
00042       if (obj->IsA() == TTree::Class()) {
00043         obj->Print();
00044       }
00045     }
00046     return;
00047   }
00048 
00049   // number of entries in a tree
00050   Long64_t numEntries(TFile *hdl, std::string const& trname) {
00051     TTree *tree = (TTree*)hdl->Get(trname.c_str());
00052     if (tree) {
00053       return tree->GetEntries();
00054     } else {
00055       std::cout << "ERR cannot find a TTree named \"" << trname << "\""
00056                 << std::endl;
00057       return -1;
00058     }
00059   }
00060 
00061   namespace {
00062     void addBranchSizes(TBranch *branch, Long64_t& size) {
00063       size += branch->GetTotalSize(); // Includes size of branch metadata
00064       // Now recurse through any subbranches.
00065       Long64_t nB = branch->GetListOfBranches()->GetEntries();
00066       for (Long64_t i = 0; i < nB; ++i) {
00067         TBranch *btemp = (TBranch *)branch->GetListOfBranches()->At(i);
00068         addBranchSizes(btemp, size);
00069       }
00070     }
00071   }
00072 
00073   void printBranchNames(TTree *tree) {
00074     if (tree != 0) {
00075       Long64_t nB = tree->GetListOfBranches()->GetEntries();
00076       for (Long64_t i = 0; i < nB; ++i) {
00077         Long64_t size = 0LL;
00078         TBranch *btemp = (TBranch *)tree->GetListOfBranches()->At(i);
00079         addBranchSizes(btemp, size);
00080         std::cout << "Branch " << i << " of " << tree->GetName() << " tree: " << btemp->GetName() << " Total size = " << size << std::endl;
00081       }
00082     } else {
00083       std::cout << "Missing Events tree?\n";
00084     }
00085   }
00086 
00087   void longBranchPrint(TTree *tr) {
00088     if (tr != 0) {
00089       Long64_t nB = tr->GetListOfBranches()->GetEntries();
00090       for (Long64_t i = 0; i < nB; ++i) {
00091         tr->GetListOfBranches()->At(i)->Print();
00092       }
00093     } else {
00094       std::cout << "Missing Events tree?\n";
00095     }
00096   }
00097 
00098   std::string getUuid(TTree *uuidTree) {
00099     FileID fid;
00100     FileID *fidPtr = &fid;
00101     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
00102     uuidTree->GetEntry(0);
00103     return fid.fid();
00104   }
00105 
00106   void printUuids(TTree *uuidTree) {
00107     FileID fid;
00108     FileID *fidPtr = &fid;
00109     uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
00110     uuidTree->GetEntry(0);
00111 
00112     std::cout << "UUID: " << fid.fid() << std::endl;
00113   }
00114 
00115   static void preIndexIntoFilePrintEventLists(TFile*, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00116     FileIndex fileIndex;
00117     FileIndex *findexPtr = &fileIndex;
00118     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != 0) {
00119       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
00120       fndx->SetAddress(&findexPtr);
00121       fndx->GetEntry(0);
00122     } else {
00123       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
00124                    "this indicates a problem with the file.  This condition should be expected with\n"
00125         "files created with earlier releases and printout of the event list will fail.\n";
00126       return;
00127     }
00128 
00129     std::cout << "\n" << fileIndex;
00130 
00131     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
00132     if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
00133     else std::cout << "This version does not support fast copy\n";
00134 
00135     if (fileIndex.allEventsInEntryOrder()) {
00136       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
00137     } else {
00138       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
00139     }
00140 
00141     fileIndex.sortBy_Run_Lumi_EventEntry();
00142     if (fileIndex.allEventsInEntryOrder()) {
00143       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
00144     } else {
00145       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
00146     }
00147     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
00148   }
00149 
00150   static void postIndexIntoFilePrintEventLists(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00151     IndexIntoFile indexIntoFile;
00152     IndexIntoFile *findexPtr = &indexIntoFile;
00153     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != 0) {
00154       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
00155       fndx->SetAddress(&findexPtr);
00156       fndx->GetEntry(0);
00157     } else {
00158       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
00159                    "this indicates a problem with the file.  This condition should be expected with\n"
00160                    "files created with earlier releases and printout of the event list will fail.\n";
00161       return;
00162     }
00163     //need to read event # from the EventAuxiliary branch
00164     TTree* eventsTree = dynamic_cast<TTree*>(tfl->Get(poolNames::eventTreeName().c_str()));
00165     TBranch* eventAuxBranch = 0;
00166     assert(0 != eventsTree);
00167     char const* const kEventAuxiliaryBranchName = "EventAuxiliary";
00168     if(eventsTree->FindBranch(kEventAuxiliaryBranchName) != 0){
00169       eventAuxBranch = eventsTree->GetBranch(kEventAuxiliaryBranchName);
00170     } else {
00171       std::cout << "Failed to find " << kEventAuxiliaryBranchName << " branch in Events TTree.  Something is wrong with this file." << std::endl;
00172       return;
00173     }
00174     EventAuxiliary eventAuxiliary;
00175     EventAuxiliary* eAPtr = &eventAuxiliary;
00176     eventAuxBranch->SetAddress(&eAPtr);
00177     std::cout << "\nPrinting IndexIntoFile contents.  This includes a list of all Runs, LuminosityBlocks\n"
00178        << "and Events stored in the root file.\n\n";
00179     std::cout << std::setw(15) << "Run"
00180        << std::setw(15) << "Lumi"
00181        << std::setw(15) << "Event"
00182        << std::setw(15) << "TTree Entry"
00183        << "\n";
00184 
00185     for(IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
00186                                         itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
00187                                         it != itEnd; ++it) {
00188       IndexIntoFile::EntryType t = it.getEntryType();
00189       std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
00190       EventNumber_t eventNum = 0;
00191       std::string type;
00192       switch(t) {
00193         case IndexIntoFile::kRun:
00194           type = "(Run)";
00195         break;
00196         case IndexIntoFile::kLumi:
00197           type = "(Lumi)";
00198         break;
00199         case IndexIntoFile::kEvent:
00200           eventAuxBranch->GetEntry(it.entry());
00201           eventNum = eventAuxiliary.id().event();
00202         break;
00203         default:
00204         break;
00205       }
00206       std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() << " " << type << std::endl;
00207     }
00208 
00209     std::cout << "\nFileFormatVersion = " << fileFormatVersion << ".  ";
00210     if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
00211     else std::cout << "This version does not support fast copy\n";
00212 
00213     if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::firstAppearanceOrder)) {
00214       std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
00215     } else {
00216       std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
00217     }
00218 
00219     // This will not work unless the other nonpersistent parts of the Index are filled first
00220     // I did not have time to implement this yet.
00221     // if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::numericalOrder)) {
00222     //   std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
00223     // } else {
00224     //   std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
00225     // }
00226     std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
00227   }
00228 
00229   void printEventLists(TFile *tfl) {
00230     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
00231     assert(0 != metaDataTree);
00232 
00233     FileFormatVersion fileFormatVersion;
00234     FileFormatVersion *fftPtr = &fileFormatVersion;
00235     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != 0) {
00236       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
00237       fft->SetAddress(&fftPtr);
00238       fft->GetEntry(0);
00239     }
00240     if(fileFormatVersion.hasIndexIntoFile()) {
00241       postIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
00242     } else {
00243       preIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
00244     }
00245   }
00246 
00247   static void preIndexIntoFilePrintEventsInLumis(TFile*, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00248     FileIndex fileIndex;
00249     FileIndex *findexPtr = &fileIndex;
00250     if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != 0) {
00251       TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
00252       fndx->SetAddress(&findexPtr);
00253       fndx->GetEntry(0);
00254     } else {
00255       std::cout << "FileIndex not found.  If this input file was created with release 1_8_0 or later\n"
00256       "this indicates a problem with the file.  This condition should be expected with\n"
00257       "files created with earlier releases and printout of the event list will fail.\n";
00258       return;
00259     }
00260 
00261     std::cout <<"\n"<< std::setw(15) << "Run"
00262     << std::setw(15) << "Lumi"
00263     << std::setw(15) << "# Events"
00264     << "\n";
00265     unsigned long nEvents = 0;
00266     unsigned long runID = 0;
00267     unsigned long lumiID = 0;
00268     for(std::vector<FileIndex::Element>::const_iterator it = fileIndex.begin(), itEnd = fileIndex.end(); it != itEnd; ++it) {
00269       if(it->getEntryType() == FileIndex::kEvent) {
00270         ++nEvents;
00271       }
00272       else if(it->getEntryType() == FileIndex::kLumi) {
00273         if(runID !=it->run_ || lumiID !=it->lumi_) {
00274           //print the previous one
00275           if(lumiID !=0) {
00276             std::cout << std::setw(15) << runID
00277             << std::setw(15) << lumiID
00278             << std::setw(15) << nEvents<<"\n";
00279           }
00280           nEvents=0;
00281           runID = it->run_;
00282           lumiID = it->lumi_;
00283         }
00284       }
00285     }
00286     //print the last one
00287     if(lumiID !=0) {
00288       std::cout << std::setw(15) << runID
00289       << std::setw(15) << lumiID
00290       << std::setw(15) << nEvents<<"\n";
00291     }
00292 
00293     std::cout << "\n";
00294   }
00295 
00296   static void postIndexIntoFilePrintEventsInLumis(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00297     IndexIntoFile indexIntoFile;
00298     IndexIntoFile *findexPtr = &indexIntoFile;
00299     if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != 0) {
00300       TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
00301       fndx->SetAddress(&findexPtr);
00302       fndx->GetEntry(0);
00303     } else {
00304       std::cout << "IndexIntoFile not found.  If this input file was created with release 1_8_0 or later\n"
00305       "this indicates a problem with the file.  This condition should be expected with\n"
00306       "files created with earlier releases and printout of the event list will fail.\n";
00307       return;
00308     }
00309     std::cout <<"\n"<< std::setw(15) << "Run"
00310     << std::setw(15) << "Lumi"
00311     << std::setw(15) << "# Events"
00312     << "\n";
00313 
00314     unsigned long nEvents = 0;
00315     unsigned long runID = 0;
00316     unsigned long lumiID = 0;
00317 
00318     for(IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
00319         itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
00320         it != itEnd; ++it) {
00321       IndexIntoFile::EntryType t = it.getEntryType();
00322       switch(t) {
00323         case IndexIntoFile::kRun:
00324           break;
00325         case IndexIntoFile::kLumi:
00326           if(runID != it.run() || lumiID != it.lumi()) {
00327             //print the previous one
00328             if(lumiID !=0) {
00329               std::cout << std::setw(15) << runID
00330               << std::setw(15) << lumiID
00331               << std::setw(15) << nEvents<<"\n";
00332             }
00333             nEvents=0;
00334             runID = it.run();
00335             lumiID = it.lumi();
00336           }
00337           break;
00338         case IndexIntoFile::kEvent:
00339           ++nEvents;
00340           break;
00341         default:
00342           break;
00343       }
00344     }
00345     //print the last one
00346     if(lumiID !=0) {
00347       std::cout << std::setw(15) << runID
00348       << std::setw(15) << lumiID
00349       << std::setw(15) << nEvents<<"\n";
00350     }
00351     std::cout << "\n";
00352   }
00353 
00354   void printEventsInLumis(TFile *tfl) {
00355     TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
00356     assert(0 != metaDataTree);
00357 
00358     FileFormatVersion fileFormatVersion;
00359     FileFormatVersion *fftPtr = &fileFormatVersion;
00360     if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != 0) {
00361       TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
00362       fft->SetAddress(&fftPtr);
00363       fft->GetEntry(0);
00364     }
00365     if(fileFormatVersion.hasIndexIntoFile()) {
00366       postIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
00367     } else {
00368       preIndexIntoFilePrintEventsInLumis(tfl, fileFormatVersion, metaDataTree);
00369     }
00370   }
00371 }