00001
00002 #include "RootOutputTree.h"
00003
00004 #include "DataFormats/Common/interface/RefCoreStreamer.h"
00005 #include "DataFormats/Provenance/interface/BranchDescription.h"
00006 #include "DataFormats/Provenance/interface/WrapperInterfaceBase.h"
00007 #include "FWCore/MessageLogger/interface/JobReport.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "FWCore/ServiceRegistry/interface/Service.h"
00010 #include "FWCore/Utilities/interface/Algorithms.h"
00011 #include "FWCore/Utilities/interface/EDMException.h"
00012 #include "FWCore/Utilities/interface/RootHandlers.h"
00013
00014 #include "TBranch.h"
00015 #include "TBranchElement.h"
00016 #include "TCollection.h"
00017 #include "TFile.h"
00018 #include "TTreeCloner.h"
00019 #include "Rtypes.h"
00020 #include "RVersion.h"
00021
00022 #include "boost/bind.hpp"
00023 #include <limits>
00024
00025 namespace edm {
00026
00027 RootOutputTree::RootOutputTree(
00028 boost::shared_ptr<TFile> filePtr,
00029 BranchType const& branchType,
00030 int splitLevel,
00031 int treeMaxVirtualSize) :
00032 filePtr_(filePtr),
00033 tree_(makeTTree(filePtr.get(), BranchTypeToProductTreeName(branchType), splitLevel)),
00034 producedBranches_(),
00035 readBranches_(),
00036 auxBranches_(),
00037 unclonedReadBranches_(),
00038 clonedReadBranchNames_(),
00039 currentlyFastCloning_(),
00040 fastCloneAuxBranches_(false) {
00041
00042 if(treeMaxVirtualSize >= 0) tree_->SetMaxVirtualSize(treeMaxVirtualSize);
00043 }
00044
00045 TTree*
00046 RootOutputTree::assignTTree(TFile* filePtr, TTree* tree) {
00047 tree->SetDirectory(filePtr);
00048
00049
00050 tree->SetAutoSave(std::numeric_limits<Long64_t>::max());
00051 return tree;
00052 }
00053
00054 TTree*
00055 RootOutputTree::makeTTree(TFile* filePtr, std::string const& name, int splitLevel) {
00056 TTree* tree = new TTree(name.c_str(), "", splitLevel);
00057 if(!tree) throw edm::Exception(errors::FatalRootError)
00058 << "Failed to create the tree: " << name << "\n";
00059 if(tree->IsZombie())
00060 throw edm::Exception(errors::FatalRootError)
00061 << "Tree: " << name << " is a zombie." << "\n";
00062
00063 return assignTTree(filePtr, tree);
00064 }
00065
00066 bool
00067 RootOutputTree::checkSplitLevelsAndBasketSizes(TTree* inputTree) const {
00068
00069 assert (inputTree != 0);
00070
00071
00072 for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end();
00073 it != itEnd; ++it) {
00074
00075 TBranch* outputBranch = *it;
00076 if(outputBranch != 0) {
00077 TBranch* inputBranch = inputTree->GetBranch(outputBranch->GetName());
00078
00079 if(inputBranch != 0) {
00080 if(inputBranch->GetSplitLevel() != outputBranch->GetSplitLevel() ||
00081 inputBranch->GetBasketSize() != outputBranch->GetBasketSize()) {
00082 return false;
00083 }
00084 }
00085 }
00086 }
00087 return true;
00088 }
00089
00090 namespace {
00091 bool checkMatchingBranches(TBranchElement* inputBranch, TBranchElement* outputBranch) {
00092 if(inputBranch->GetStreamerType() != outputBranch->GetStreamerType()) {
00093 return false;
00094 }
00095 TObjArray* inputArray = inputBranch->GetListOfBranches();
00096 TObjArray* outputArray = outputBranch->GetListOfBranches();
00097
00098 if(outputArray->GetSize() < inputArray->GetSize()) {
00099 return false;
00100 }
00101 TIter iter(outputArray);
00102 TObject* obj = 0;
00103 while((obj = iter.Next()) != 0) {
00104 TBranchElement* outBranch = dynamic_cast<TBranchElement*>(obj);
00105 if(outBranch) {
00106 TBranchElement* inBranch = dynamic_cast<TBranchElement*>(inputArray->FindObject(outBranch->GetName()));
00107 if(!inBranch) {
00108 return false;
00109 }
00110 if(!checkMatchingBranches(inBranch, outBranch)) {
00111 return false;
00112 }
00113 }
00114 }
00115 return true;
00116 }
00117 }
00118
00119 bool RootOutputTree::checkIfFastClonable(TTree* inputTree) const {
00120
00121 if(inputTree == 0) return false;
00122
00123
00124 for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end(); it != itEnd; ++it) {
00125 TBranchElement* outputBranch = dynamic_cast<TBranchElement*>(*it);
00126 if(outputBranch != 0) {
00127 TBranchElement* inputBranch = dynamic_cast<TBranchElement*>(inputTree->GetBranch(outputBranch->GetName()));
00128 if(inputBranch != 0) {
00129
00130 if(!checkMatchingBranches(inputBranch, outputBranch)) {
00131 LogInfo("FastCloning")
00132 << "Fast Cloning disabled because a data member has been added to split branch: " << inputBranch->GetName() << "\n.";
00133 return false;
00134 }
00135 }
00136 }
00137 }
00138 return true;
00139 }
00140
00141 bool RootOutputTree::checkEntriesInReadBranches(Long64_t expectedNumberOfEntries) const {
00142 for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end(); it != itEnd; ++it) {
00143 if((*it)->GetEntries() != expectedNumberOfEntries) {
00144 return false;
00145 }
00146 }
00147 return true;
00148 }
00149
00150 void
00151 RootOutputTree::fastCloneTTree(TTree* in, std::string const& option) {
00152 if(in->GetEntries() != 0) {
00153 TObjArray* branches = tree_->GetListOfBranches();
00154
00155
00156 std::map<Int_t, TBranch *> auxIndexes;
00157 bool mustRemoveSomeAuxs = false;
00158 if(!fastCloneAuxBranches_) {
00159 for(std::vector<TBranch *>::const_iterator it = auxBranches_.begin(), itEnd = auxBranches_.end();
00160 it != itEnd; ++it) {
00161 int auxIndex = branches->IndexOf(*it);
00162 assert (auxIndex >= 0);
00163 auxIndexes.insert(std::make_pair(auxIndex, *it));
00164 branches->RemoveAt(auxIndex);
00165 }
00166 mustRemoveSomeAuxs = true;
00167 }
00168
00169
00170 for(std::vector<TBranch *>::const_iterator it = unclonedAuxBranches_.begin(),
00171 itEnd = unclonedAuxBranches_.end();
00172 it != itEnd; ++it) {
00173 int auxIndex = branches->IndexOf(*it);
00174 assert (auxIndex >= 0);
00175 auxIndexes.insert(std::make_pair(auxIndex, *it));
00176 branches->RemoveAt(auxIndex);
00177 mustRemoveSomeAuxs = true;
00178 }
00179
00180 if(mustRemoveSomeAuxs) {
00181 branches->Compress();
00182 }
00183
00184 TTreeCloner cloner(in, tree_, option.c_str(), TTreeCloner::kNoWarnings|TTreeCloner::kIgnoreMissingTopLevel);
00185
00186 if(!cloner.IsValid()) {
00187
00188 static const char* okerror = "One of the export branch";
00189 if(strncmp(cloner.GetWarning(), okerror, strlen(okerror)) == 0) {
00190
00191 } else {
00192 throw edm::Exception(errors::FatalRootError)
00193 << "invalid TTreeCloner (" << cloner.GetWarning() << ")\n";
00194 }
00195 }
00196 tree_->SetEntries(tree_->GetEntries() + in->GetEntries());
00197 Service<RootHandlers> rootHandler;
00198 rootHandler->enableErrorHandlerWithoutWarnings();
00199 cloner.Exec();
00200 rootHandler->enableErrorHandler();
00201 if(mustRemoveSomeAuxs) {
00202 for(std::map<Int_t, TBranch *>::const_iterator it = auxIndexes.begin(), itEnd = auxIndexes.end();
00203 it != itEnd; ++it) {
00204
00205 Int_t last = branches->GetLast();
00206 if(last >= 0) {
00207 branches->AddAtAndExpand(branches->At(last), last+1);
00208 for(Int_t ind = last-1; ind >= it->first; --ind) {
00209 branches->AddAt(branches->At(ind), ind+1);
00210 };
00211 branches->AddAt(it->second, it->first);
00212 } else {
00213 branches->Add(it->second);
00214 }
00215 }
00216 }
00217 }
00218 }
00219
00220 void
00221 RootOutputTree::writeTTree(TTree* tree) {
00222 if(tree->GetNbranches() != 0) {
00223 tree->SetEntries(-1);
00224 }
00225 setRefCoreStreamer(true);
00226 tree->AutoSave("FlushBaskets");
00227 }
00228
00229 void
00230 RootOutputTree::fillTTree(std::vector<TBranch*> const& branches) {
00231 for_all(branches, boost::bind(&TBranch::Fill, _1));
00232 }
00233
00234 void
00235 RootOutputTree::writeTree() const {
00236 writeTTree(tree_);
00237 }
00238
00239 void
00240 RootOutputTree::maybeFastCloneTree(bool canFastClone, bool canFastCloneAux, TTree* tree, std::string const& option) {
00241 unclonedReadBranches_.clear();
00242 clonedReadBranchNames_.clear();
00243 currentlyFastCloning_ = canFastClone && !readBranches_.empty();
00244 if(currentlyFastCloning_) {
00245 fastCloneAuxBranches_ = canFastCloneAux;
00246 fastCloneTTree(tree, option);
00247 for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end();
00248 it != itEnd; ++it) {
00249 if((*it)->GetEntries() == tree_->GetEntries()) {
00250 clonedReadBranchNames_.insert(std::string((*it)->GetName()));
00251 } else {
00252 unclonedReadBranches_.push_back(*it);
00253 }
00254 }
00255 Service<JobReport> reportSvc;
00256 reportSvc->reportFastClonedBranches(clonedReadBranchNames_, tree_->GetEntries());
00257 }
00258 }
00259
00260 void
00261 RootOutputTree::fillTree() const {
00262 if(currentlyFastCloning_) {
00263 if(!fastCloneAuxBranches_)fillTTree(auxBranches_);
00264 fillTTree(unclonedAuxBranches_);
00265 fillTTree(producedBranches_);
00266 fillTTree(unclonedReadBranches_);
00267 } else {
00268 tree_->Fill();
00269 }
00270 }
00271
00272 void
00273 RootOutputTree::addBranch(std::string const& branchName,
00274 std::string const& className,
00275 WrapperInterfaceBase const* interface,
00276 void const*& pProd,
00277 int splitLevel,
00278 int basketSize,
00279 bool produced) {
00280 assert(splitLevel != BranchDescription::invalidSplitLevel);
00281 assert(basketSize != BranchDescription::invalidBasketSize);
00282 TBranch* branch = tree_->Branch(branchName.c_str(),
00283 className.c_str(),
00284 &pProd,
00285 basketSize,
00286 splitLevel);
00287 assert(branch != 0);
00288 if(pProd != 0) {
00289
00290 interface->deleteProduct(pProd);
00291 pProd = 0;
00292 }
00293 if(produced) {
00294 producedBranches_.push_back(branch);
00295 } else {
00296 readBranches_.push_back(branch);
00297 }
00298 }
00299
00300 void
00301 RootOutputTree::close() {
00302
00303
00304 auxBranches_.clear();
00305 unclonedAuxBranches_.clear();
00306 producedBranches_.clear();
00307 readBranches_.clear();
00308 unclonedReadBranches_.clear();
00309 tree_ = 0;
00310 filePtr_.reset();
00311 }
00312 }