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/FileID.h"
00005 #include "DataFormats/Provenance/interface/FileFormatVersion.h"
00006 #include "DataFormats/Provenance/interface/FileIndex.h"
00007 #include "DataFormats/Provenance/interface/IndexIntoFile.h"
00008 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
00009
00010 #include <iostream>
00011 #include <iomanip>
00012
00013 #include "TFile.h"
00014 #include "TList.h"
00015 #include "TIterator.h"
00016 #include "TKey.h"
00017 #include "TTree.h"
00018 #include "TObject.h"
00019 #include "TBranch.h"
00020
00021 namespace edm {
00022
00023
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
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
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 void printBranchNames(TTree *tree) {
00062 if (tree != 0) {
00063 Long64_t nB = tree->GetListOfBranches()->GetEntries();
00064 for (Long64_t i = 0; i < nB; ++i) {
00065 TBranch *btemp = (TBranch *)tree->GetListOfBranches()->At(i);
00066 std::cout << "Branch " << i <<" of " << tree->GetName() <<" tree: " << btemp->GetName() << " Total size = " << btemp->GetTotalSize() << std::endl;
00067 }
00068 } else {
00069 std::cout << "Missing Events tree?\n";
00070 }
00071 }
00072
00073 void longBranchPrint(TTree *tr) {
00074 if (tr != 0) {
00075 Long64_t nB = tr->GetListOfBranches()->GetEntries();
00076 for (Long64_t i = 0; i < nB; ++i) {
00077 tr->GetListOfBranches()->At(i)->Print();
00078 }
00079 } else {
00080 std::cout << "Missing Events tree?\n";
00081 }
00082 }
00083
00084 std::string getUuid(TTree *uuidTree) {
00085 FileID fid;
00086 FileID *fidPtr = &fid;
00087 uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
00088 uuidTree->GetEntry(0);
00089 return fid.fid();
00090 }
00091
00092 void printUuids(TTree *uuidTree) {
00093 FileID fid;
00094 FileID *fidPtr = &fid;
00095 uuidTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
00096 uuidTree->GetEntry(0);
00097
00098 std::cout << "UUID: " << fid.fid() << std::endl;
00099 }
00100
00101 static void preIndexIntoFilePrintEventLists(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00102 FileIndex fileIndex;
00103 FileIndex *findexPtr = &fileIndex;
00104 if (metaDataTree->FindBranch(poolNames::fileIndexBranchName().c_str()) != 0) {
00105 TBranch *fndx = metaDataTree->GetBranch(poolNames::fileIndexBranchName().c_str());
00106 fndx->SetAddress(&findexPtr);
00107 fndx->GetEntry(0);
00108 } else {
00109 std::cout << "FileIndex not found. If this input file was created with release 1_8_0 or later\n"
00110 "this indicates a problem with the file. This condition should be expected with\n"
00111 "files created with earlier releases and printout of the event list will fail.\n";
00112 return;
00113 }
00114
00115 std::cout << "\n" << fileIndex;
00116
00117 std::cout << "\nFileFormatVersion = " << fileFormatVersion << ". ";
00118 if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
00119 else std::cout << "This version does not support fast copy\n";
00120
00121 if (fileIndex.allEventsInEntryOrder()) {
00122 std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = False\" mode\n";
00123 } else {
00124 std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = False\" mode\n";
00125 }
00126
00127 fileIndex.sortBy_Run_Lumi_EventEntry();
00128 if (fileIndex.allEventsInEntryOrder()) {
00129 std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort\" mode\n";
00130 } else {
00131 std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort\" mode\n";
00132 }
00133 std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
00134 }
00135
00136 static void postIndexIntoFilePrintEventLists(TFile* tfl, FileFormatVersion const& fileFormatVersion, TTree *metaDataTree) {
00137 IndexIntoFile indexIntoFile;
00138 IndexIntoFile *findexPtr = &indexIntoFile;
00139 if (metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != 0) {
00140 TBranch *fndx = metaDataTree->GetBranch(poolNames::indexIntoFileBranchName().c_str());
00141 fndx->SetAddress(&findexPtr);
00142 fndx->GetEntry(0);
00143 } else {
00144 std::cout << "IndexIntoFile not found. If this input file was created with release 1_8_0 or later\n"
00145 "this indicates a problem with the file. This condition should be expected with\n"
00146 "files created with earlier releases and printout of the event list will fail.\n";
00147 return;
00148 }
00149
00150 TTree* eventsTree = dynamic_cast<TTree*>(tfl->Get(edm::poolNames::eventTreeName().c_str()));
00151 TBranch* eventAuxBranch = 0;
00152 assert(0 != eventsTree);
00153 char const* const kEventAuxiliaryBranchName = "EventAuxiliary";
00154 if(eventsTree->FindBranch(kEventAuxiliaryBranchName) != 0){
00155 eventAuxBranch = eventsTree->GetBranch(kEventAuxiliaryBranchName);
00156 } else {
00157 std::cout << "Failed to find " << kEventAuxiliaryBranchName << " branch in Events TTree. Something is wrong with this file." << std::endl;
00158 return;
00159 }
00160 EventAuxiliary eventAuxiliary;
00161 EventAuxiliary* eAPtr = &eventAuxiliary;
00162 eventAuxBranch->SetAddress(&eAPtr);
00163 std::cout << "\nPrinting IndexIntoFile contents. This includes a list of all Runs, LuminosityBlocks\n"
00164 << "and Events stored in the root file.\n\n";
00165 std::cout << std::setw(15) << "Run"
00166 << std::setw(15) << "Lumi"
00167 << std::setw(15) << "Event"
00168 << std::setw(15) << "TTree Entry"
00169 << "\n";
00170
00171 for(IndexIntoFile::IndexIntoFileItr it = indexIntoFile.begin(IndexIntoFile::firstAppearanceOrder),
00172 itEnd = indexIntoFile.end(IndexIntoFile::firstAppearanceOrder);
00173 it != itEnd; ++it) {
00174 IndexIntoFile::EntryType t = it.getEntryType();
00175 std::cout << std::setw(15) << it.run() << std::setw(15) << it.lumi();
00176 EventNumber_t eventNum = 0;
00177 std::string type;
00178 switch(t) {
00179 case IndexIntoFile::kRun:
00180 type = "(Run)";
00181 break;
00182 case IndexIntoFile::kLumi:
00183 type = "(Lumi)";
00184 break;
00185 case IndexIntoFile::kEvent:
00186 eventAuxBranch->GetEntry(it.entry());
00187 eventNum = eventAuxiliary.id().event();
00188 break;
00189 default:
00190 break;
00191 }
00192 std::cout << std::setw(15) << eventNum << std::setw(15) << it.entry() << " " << type << std::endl;
00193 }
00194
00195 std::cout << "\nFileFormatVersion = " << fileFormatVersion << ". ";
00196 if (fileFormatVersion.fastCopyPossible()) std::cout << "This version supports fast copy\n";
00197 else std::cout << "This version does not support fast copy\n";
00198
00199 if (indexIntoFile.iterationWillBeInEntryOrder(IndexIntoFile::firstAppearanceOrder)) {
00200 std::cout << "Events are sorted such that fast copy is possible in the \"noEventSort = false\" mode\n";
00201 } else {
00202 std::cout << "Events are sorted such that fast copy is NOT possible in the \"noEventSort = false\" mode\n";
00203 }
00204
00205
00206
00207
00208
00209
00210
00211
00212 std::cout << "(Note that other factors can prevent fast copy from occurring)\n\n";
00213 }
00214
00215 void printEventLists(TFile *tfl) {
00216 TTree *metaDataTree = dynamic_cast<TTree *>(tfl->Get(poolNames::metaDataTreeName().c_str()));
00217
00218 FileFormatVersion fileFormatVersion;
00219 FileFormatVersion *fftPtr = &fileFormatVersion;
00220 if (metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != 0) {
00221 TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
00222 fft->SetAddress(&fftPtr);
00223 fft->GetEntry(0);
00224 }
00225 if(fileFormatVersion.hasIndexIntoFile()) {
00226 postIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
00227 } else {
00228 preIndexIntoFilePrintEventLists(tfl, fileFormatVersion, metaDataTree);
00229 }
00230 }
00231 }