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