CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RootOutputTree.cc
Go to the documentation of this file.
1 
2 #include "RootOutputTree.h"
3 
12 
13 #include "TBranch.h"
14 #include "TBranchElement.h"
15 #include "TCollection.h"
16 #include "TFile.h"
17 #include "TTreeCloner.h"
18 #include "Rtypes.h"
19 #include "RVersion.h"
20 
21 #include "boost/bind.hpp"
22 #include <limits>
23 
24 namespace edm {
25 
27  boost::shared_ptr<TFile> filePtr,
28  BranchType const& branchType,
29  int bufSize,
30  int splitLevel,
31  int treeMaxVirtualSize) :
32  filePtr_(filePtr),
33  tree_(makeTTree(filePtr.get(), BranchTypeToProductTreeName(branchType), splitLevel)),
34  producedBranches_(),
35  readBranches_(),
36  auxBranches_(),
37  unclonedReadBranches_(),
38  unclonedReadBranchNames_(),
39  currentlyFastCloning_(),
40  fastCloneAuxBranches_(false) {
41 
42  if(treeMaxVirtualSize >= 0) tree_->SetMaxVirtualSize(treeMaxVirtualSize);
43  }
44 
45  TTree*
46  RootOutputTree::assignTTree(TFile* filePtr, TTree* tree) {
47  tree->SetDirectory(filePtr);
48  // Turn off autosaving because it is such a memory hog and we are not using
49  // this check-pointing feature anyway.
50  tree->SetAutoSave(std::numeric_limits<Long64_t>::max());
51  return tree;
52  }
53 
54  TTree*
55  RootOutputTree::makeTTree(TFile* filePtr, std::string const& name, int splitLevel) {
56  TTree* tree = new TTree(name.c_str(), "", splitLevel);
57  if(!tree) throw edm::Exception(errors::FatalRootError)
58  << "Failed to create the tree: " << name << "\n";
59  if(tree->IsZombie())
61  << "Tree: " << name << " is a zombie." << "\n";
62 
63  return assignTTree(filePtr, tree);
64  }
65 
66  bool
68 
69  assert (inputTree != 0);
70 
71  // Do the split level and basket size match in the input and output?
72  for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end();
73  it != itEnd; ++it) {
74 
75  TBranch* outputBranch = *it;
76  if(outputBranch != 0) {
77  TBranch* inputBranch = inputTree->GetBranch(outputBranch->GetName());
78 
79  if(inputBranch != 0) {
80  if(inputBranch->GetSplitLevel() != outputBranch->GetSplitLevel() ||
81  inputBranch->GetBasketSize() != outputBranch->GetBasketSize()) {
82  return false;
83  }
84  }
85  }
86  }
87  return true;
88  }
89 
90  namespace {
91  bool checkMatchingBranches(TBranchElement* inputBranch, TBranchElement* outputBranch) {
92  if(inputBranch->GetStreamerType() != outputBranch->GetStreamerType()) {
93  return false;
94  }
95  TObjArray* inputArray = inputBranch->GetListOfBranches();
96  TObjArray* outputArray = outputBranch->GetListOfBranches();
97 
98  if(outputArray->GetSize() < inputArray->GetSize()) {
99  return false;
100  }
101  TIter iter(outputArray);
102  TObject* obj = 0;
103  while((obj = iter.Next()) != 0) {
104  TBranchElement* outBranch = dynamic_cast<TBranchElement*>(obj);
105  if(outBranch) {
106  TBranchElement* inBranch = dynamic_cast<TBranchElement*>(inputArray->FindObject(outBranch->GetName()));
107  if(!inBranch) {
108  return false;
109  }
110  if(!checkMatchingBranches(inBranch, outBranch)) {
111  return false;
112  }
113  }
114  }
115  return true;
116  }
117  }
118 
119  bool RootOutputTree::checkIfFastClonable(TTree* inputTree) const {
120 
121  if(inputTree == 0) return false;
122 
123  // Do the sub-branches match in the input and output. Extra sub-branches in the input are OK for fast cloning, but not in the output.
124  for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end(); it != itEnd; ++it) {
125  TBranchElement* outputBranch = dynamic_cast<TBranchElement*>(*it);
126  if(outputBranch != 0) {
127  TBranchElement* inputBranch = dynamic_cast<TBranchElement*>(inputTree->GetBranch(outputBranch->GetName()));
128  if(inputBranch != 0) {
129  // We have a matching top level branch. Do the recursive check on subbranches.
130  if(!checkMatchingBranches(inputBranch, outputBranch)) {
131  LogInfo("FastCloning")
132  << "Fast Cloning disabled because a data member has been added to split branch: " << inputBranch->GetName() << "\n.";
133  }
134  }
135  }
136  }
137  return true;
138  }
139 
140  void
141  RootOutputTree::fastCloneTTree(TTree* in, std::string const& option) {
142  if(in->GetEntries() != 0) {
143  TObjArray* branches = tree_->GetListOfBranches();
144  // If any products were produced (not just event products), the EventAuxiliary will be modified.
145  // In that case, don't fast copy auxiliary branches. Remove them, and add back after fast copying.
146  std::map<Int_t, TBranch *> auxIndexes;
147  if (!fastCloneAuxBranches_) {
148  for (std::vector<TBranch *>::const_iterator it = auxBranches_.begin(), itEnd = auxBranches_.end();
149  it != itEnd; ++it) {
150  int auxIndex = branches->IndexOf(*it);
151  assert (auxIndex >= 0);
152  auxIndexes.insert(std::make_pair(auxIndex, *it));
153  branches->RemoveAt(auxIndex);
154  }
155  branches->Compress();
156  }
157 
158 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
159  TTreeCloner cloner(in, tree_, option.c_str(), TTreeCloner::kNoWarnings|TTreeCloner::kIgnoreMissingTopLevel);
160 #else
161  TTreeCloner cloner(in, tree_, option.c_str());
162 #endif
163 
164  if(!cloner.IsValid()) {
165 
166 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
167  // Let's check why
168  static const char* okerror = "One of the export branch";
169  if ( strncmp(cloner.GetWarning(),okerror,strlen(okerror)) == 0 ) {
170  // That's fine we will handle it;
171  }
172  else {
174  << "invalid TTreeCloner\n";
175  }
176 #else
178  << "invalid TTreeCloner\n";
179 #endif
180  }
181  tree_->SetEntries(tree_->GetEntries() + in->GetEntries());
182  Service<RootHandlers> rootHandler;
183  rootHandler->enableErrorHandlerWithoutWarnings();
184  cloner.Exec();
185  rootHandler->enableErrorHandler();
186  if (!fastCloneAuxBranches_) {
187  for (std::map<Int_t, TBranch *>::const_iterator it = auxIndexes.begin(), itEnd = auxIndexes.end();
188  it != itEnd; ++it) {
189  // Add the auxiliary branches back after fast copying the rest of the tree.
190  Int_t last = branches->GetLast();
191  if (last >= 0) {
192  branches->AddAtAndExpand(branches->At(last), last+1);
193  for(Int_t ind = last-1; ind >= it->first; --ind) {
194  branches->AddAt(branches->At(ind), ind+1);
195  };
196  branches->AddAt(it->second, it->first);
197  } else {
198  branches->Add(it->second);
199  }
200  }
201  }
202  }
203  }
204 
205  void
207  if(tree->GetNbranches() != 0) {
208  tree->SetEntries(-1);
209  }
210  setRefCoreStreamer(true);
211  tree->AutoSave("FlushBaskets");
212  }
213 
214  void
215  RootOutputTree::fillTTree(TTree* tree, std::vector<TBranch*> const& branches) {
216  for_all(branches, boost::bind(&TBranch::Fill, _1));
217  }
218 
219  void
221  writeTTree(tree_);
222  }
223 
224  void
225  RootOutputTree::maybeFastCloneTree(bool canFastClone, bool canFastCloneAux, TTree* tree, std::string const& option) {
226  unclonedReadBranches_.clear();
227  unclonedReadBranchNames_.clear();
228  currentlyFastCloning_ = canFastClone && !readBranches_.empty();
230  fastCloneAuxBranches_ = canFastCloneAux;
231  fastCloneTTree(tree, option);
232  for(std::vector<TBranch*>::const_iterator it = readBranches_.begin(), itEnd = readBranches_.end();
233  it != itEnd; ++it) {
234  if((*it)->GetEntries() != tree_->GetEntries()) {
235  unclonedReadBranches_.push_back(*it);
236  unclonedReadBranchNames_.insert(std::string((*it)->GetName()));
237  }
238  }
239  }
240  }
241 
242  void
244  if (currentlyFastCloning_) {
248  } else {
249  tree_->Fill();
250  }
251  }
252 
253  void
254  RootOutputTree::addBranch(std::string const& branchName,
255  std::string const& className,
256  void const*& pProd,
257  int splitLevel,
258  int basketSize,
259  bool produced) {
260  assert(splitLevel != BranchDescription::invalidSplitLevel);
261  assert(basketSize != BranchDescription::invalidBasketSize);
262  TBranch* branch = tree_->Branch(branchName.c_str(),
263  className.c_str(),
264  &pProd,
265  basketSize,
266  splitLevel);
267  if(produced) {
268  producedBranches_.push_back(branch);
269  } else {
270  readBranches_.push_back(branch);
271  }
272  }
273 
274  void
276  // The TFile was just closed.
277  // Just to play it safe, zero all pointers to quantities in the file.
278  auxBranches_.clear();
279  producedBranches_.clear();
280  readBranches_.clear();
281  unclonedReadBranches_.clear();
282  tree_ = 0;
283  filePtr_.reset();
284  }
285 }
void writeTree() const
TTree *const tree() const
std::set< std::string > unclonedReadBranchNames_
static int const invalidSplitLevel
void fillTree() const
static int const invalidBasketSize
bool checkSplitLevelsAndBasketSizes(TTree *inputTree) const
void setRefCoreStreamer(bool resetAll=false)
bool checkIfFastClonable(TTree *inputTree) const
BranchType
Definition: BranchType.h:11
std::vector< TBranch * > producedBranches_
std::vector< TBranch * > auxBranches_
Func for_all(ForwardSequence &s, Func f)
wrapper for std::for_each
Definition: Algorithms.h:16
boost::shared_ptr< TFile > filePtr_
RootOutputTree(boost::shared_ptr< TFile > filePtr, BranchType const &branchType, int bufSize, int splitLevel, int treeMaxVirtualSize)
static TTree * assignTTree(TFile *file, TTree *tree)
tuple obj
Example code starts here #.
Definition: VarParsing.py:655
void addBranch(std::string const &branchName, std::string const &className, void const *&pProd, int splitLevel, int basketSize, bool produced)
const T & max(const T &a, const T &b)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
static void fillTTree(TTree *tree, std::vector< TBranch * > const &branches)
std::vector< TBranch * > unclonedReadBranches_
std::string const & BranchTypeToProductTreeName(BranchType const &branchType)
Definition: BranchType.cc:96
static TTree * makeTTree(TFile *filePtr, std::string const &name, int splitLevel)
eventsetup::produce::Produce produced
Definition: ESProducts.cc:21
static void writeTTree(TTree *tree)
std::vector< TBranch * > readBranches_
void maybeFastCloneTree(bool canFastClone, bool canFastCloneAux, TTree *tree, std::string const &option)
T get(const Candidate &c)
Definition: component.h:56
std::string className(const T &t)
Definition: ClassName.h:30
void fastCloneTTree(TTree *in, std::string const &option)