CMS 3D CMS Logo

TreeReader.cc
Go to the documentation of this file.
1 #include <cstdint>
2 #include <utility>
3 #include <cstring>
4 #include <string>
5 #include <vector>
6 #include <map>
7 
8 #include <TString.h>
9 #include <TTree.h>
10 #include <TBranch.h>
11 #include <TLeaf.h>
12 #include <TList.h>
13 #include <TKey.h>
14 
16 
20 
21 namespace PhysicsTools {
22 
23  const double TreeReader::kOptVal = -999.0;
24 
26 
27  TreeReader::TreeReader(const TreeReader &orig) { this->operator=(orig); }
28 
29  TreeReader::TreeReader(TTree *tree, bool skipTarget, bool skipWeight) : tree(tree), upToDate(false) {
30  automaticAdd(skipTarget, skipWeight);
31  }
32 
34 
36  reset();
37 
38  tree = orig.tree;
39 
40  multiDouble.resize(orig.multiDouble.size());
41  multiFloat.resize(orig.multiFloat.size());
42  multiInt.resize(orig.multiInt.size());
43  multiBool.resize(orig.multiBool.size());
44 
45  singleDouble.resize(orig.singleDouble.size());
46  singleFloat.resize(orig.singleFloat.size());
47  singleInt.resize(orig.singleInt.size());
48  singleBool.resize(orig.singleBool.size());
49 
50  valueMap = orig.valueMap;
51 
52  return *this;
53  }
54 
55  void TreeReader::setTree(TTree *tree) {
56  this->tree = tree;
57  upToDate = false;
58  }
59 
60  void TreeReader::addBranch(const std::string &expression, AtomicId name, bool opt) {
61  if (!tree)
62  throw cms::Exception("NoTreeAvailable") << "No TTree set in TreeReader::addBranch." << std::endl;
63 
64  TBranch *branch = tree->GetBranch(expression.c_str());
65  if (!branch)
66  throw cms::Exception("BranchMissing") << "Tree branch \"" << expression << "\" missing." << std::endl;
67 
68  addBranch(branch, name, opt);
69  }
70 
71  void TreeReader::addBranch(TBranch *branch, AtomicId name, bool opt) {
72  TString branchName = branch->GetName();
73  if (!name)
74  name = (const char *)branchName;
75 
76  TLeaf *leaf = dynamic_cast<TLeaf *>(branch->GetLeaf(branchName));
77  if (!leaf)
78  throw cms::Exception("InvalidBranch") << "Tree branch \"" << branchName << "\" has no leaf." << std::endl;
79 
80  TString typeName = leaf->GetTypeName();
81  char typeId = 0;
82  bool multi = false;
83  if (typeName == "Double_t" || typeName == "double")
84  typeId = 'D';
85  else if (typeName == "Float_t" || typeName == "float")
86  typeId = 'F';
87  else if (typeName == "Int_t" || typeName == "int")
88  typeId = 'I';
89  else if (typeName == "Bool_t" || typeName == "bool")
90  typeId = 'B';
91  else {
92  multi = true;
93  if (typeName == "vector<double>" || typeName == "Vector<Double_t>")
94  typeId = 'D';
95  else if (typeName == "vector<float>" || typeName == "Vector<Float_t>")
96  typeId = 'F';
97  else if (typeName == "vector<int>" || typeName == "Vector<Int_t>")
98  typeId = 'I';
99  else if (typeName == "vector<bool>" || typeName == "Vector<Bool_t>")
100  typeId = 'B';
101  }
102 
103  if (!typeId)
104  throw cms::Exception("InvalidBranch") << "Tree branch \"" << branchName
105  << "\" is of "
106  "unsupported type \""
107  << typeName << "\"." << std::endl;
108 
109  if (multi)
110  addTypeMulti(name, nullptr, typeId);
111  else
112  addTypeSingle(name, nullptr, typeId, opt);
113 
114  valueMap[name].setBranchName(branch->GetName());
115  }
116 
117  void TreeReader::setOptional(AtomicId name, bool opt, double optVal) {
118  std::map<AtomicId, Value>::iterator pos = valueMap.find(name);
119  if (pos == valueMap.end())
120  throw cms::Exception("UnknownVariable") << "Variable \"" << name
121  << "\" is not known to the "
122  "TreeReader."
123  << std::endl;
124 
125  pos->second.setOpt(opt, optVal);
126  }
127 
128  void TreeReader::addTypeSingle(AtomicId name, const void *value, char type, bool opt) {
129  std::map<AtomicId, Value>::const_iterator pos = valueMap.find(name);
130  if (pos != valueMap.end())
131  throw cms::Exception("DuplicateVariable") << "Duplicate Variable \"" << name << "\"." << std::endl;
132 
133  if (type != 'D' && type != 'F' && type != 'I' && type != 'B')
134  throw cms::Exception("InvalidType") << "Unsupported type '" << type
135  << "' in call to"
136  "TreeReader::addTypeSingle."
137  << std::endl;
138 
139  int index = -1;
140  if (!value) {
141  switch (type) {
142  case 'D':
143  index = (int)singleDouble.size();
144  singleDouble.push_back(Double_t());
145  break;
146  case 'F':
147  index = (int)singleFloat.size();
148  singleFloat.push_back(Float_t());
149  break;
150  case 'I':
151  index = (int)singleInt.size();
152  singleInt.push_back(Int_t());
153  break;
154  case 'B':
155  index = (int)singleBool.size();
156  singleBool.push_back(Bool());
157  break;
158  }
159  }
160 
161  valueMap[name] = Value(index, false, opt, type);
162  if (value)
163  valueMap[name].setPtr(value);
164 
165  upToDate = false;
166  }
167 
168  template <typename T>
169  static std::pair<void *, std::vector<T> > makeMulti() {
170  return std::pair<void *, std::vector<T> >(nullptr, std::vector<T>());
171  }
172 
173  void TreeReader::addTypeMulti(AtomicId name, const void *value, char type) {
174  std::map<AtomicId, Value>::const_iterator pos = valueMap.find(name);
175  if (pos != valueMap.end())
176  throw cms::Exception("DuplicateVariable") << "Duplicate Variable \"" << name << "\"." << std::endl;
177 
178  if (type != 'D' && type != 'F' && type != 'I' && type != 'B')
179  throw cms::Exception("InvalidType") << "Unsupported type '" << type
180  << "' in call to"
181  "TreeReader::addTypeMulti."
182  << std::endl;
183 
184  int index = -1;
185  if (!value) {
186  switch (type) {
187  case 'D':
188  index = (int)multiDouble.size();
189  multiDouble.push_back(makeMulti<Double_t>());
190  break;
191  case 'F':
192  index = (int)multiFloat.size();
193  multiFloat.push_back(makeMulti<Float_t>());
194  break;
195  case 'I':
196  index = (int)multiInt.size();
197  multiInt.push_back(makeMulti<Int_t>());
198  break;
199  case 'B':
200  index = (int)multiBool.size();
201  multiBool.push_back(makeMulti<Bool_t>());
202  break;
203  }
204  }
205 
206  valueMap[name] = Value(index, true, false, type);
207  if (value)
208  valueMap[name].setPtr(value);
209 
210  upToDate = false;
211  }
212 
213  void TreeReader::automaticAdd(bool skipTarget, bool skipWeight) {
214  if (!tree)
215  throw cms::Exception("NoTreeAvailable") << "No TTree set in TreeReader::automaticAdd." << std::endl;
216 
217  TIter iter(tree->GetListOfBranches());
218  TObject *obj;
219  while ((obj = iter())) {
220  TBranch *branch = dynamic_cast<TBranch *>(obj);
221  if (!branch)
222  continue;
223 
224  if (skipTarget && !std::strcmp(branch->GetName(), "__TARGET__"))
225  continue;
226 
227  if (skipWeight && !std::strcmp(branch->GetName(), "__WEIGHT__"))
228  continue;
229 
230  addBranch(branch);
231  }
232  }
233 
235  multiDouble.clear();
236  multiFloat.clear();
237  multiInt.clear();
238  multiBool.clear();
239 
240  singleDouble.clear();
241  singleFloat.clear();
242  singleInt.clear();
243  singleBool.clear();
244 
245  valueMap.clear();
246 
247  upToDate = false;
248  }
249 
251  if (!tree)
252  throw cms::Exception("NoTreeAvailable") << "No TTree set in TreeReader::automaticAdd." << std::endl;
253 
254  for (std::map<AtomicId, Value>::iterator iter = valueMap.begin(); iter != valueMap.end(); iter++)
255  iter->second.update(this);
256 
257  upToDate = true;
258  }
259 
261  if (!tree)
262  throw cms::Exception("NoTreeAvailable") << "No TTree set in TreeReader::automaticAdd." << std::endl;
263 
264  if (!upToDate)
265  update();
266 
267  Long64_t entries = tree->GetEntries();
268  for (Long64_t entry = 0; entry < entries; entry++) {
269  tree->GetEntry(entry);
270  fill(mva);
271  }
272 
273  return entries;
274  }
275 
277  for (std::map<AtomicId, Value>::const_iterator iter = valueMap.begin(); iter != valueMap.end(); iter++)
278  iter->second.fill(iter->first, this);
279 
280  double result = mva->eval(values);
281  values.clear();
282 
283  return result;
284  }
285 
287  for (std::map<AtomicId, Value>::const_iterator iter = valueMap.begin(); iter != valueMap.end(); iter++)
288  iter->second.fill(iter->first, this);
289 
291  values.clear();
292 
293  return result;
294  }
295 
296  std::vector<AtomicId> TreeReader::variables() const {
297  std::vector<AtomicId> result;
298  for (std::map<AtomicId, Value>::const_iterator iter = valueMap.begin(); iter != valueMap.end(); iter++)
299  result.push_back(iter->first);
300 
301  return result;
302  }
303 
305  if (ptr)
306  return;
307 
308  void *value = nullptr;
309  if (multiple) {
310  switch (type) {
311  case 'D':
312  reader->multiDouble[index].first = &reader->multiDouble[index].second;
313  value = &reader->multiDouble[index].first;
314  break;
315  case 'F':
316  reader->multiFloat[index].first = &reader->multiFloat[index].second;
317  value = &reader->multiFloat[index].first;
318  break;
319  case 'I':
320  reader->multiInt[index].first = &reader->multiInt[index].second;
321  value = &reader->multiInt[index].first;
322  break;
323  case 'B':
324  reader->multiBool[index].first = value;
325  value = &reader->multiBool[index].first;
326  break;
327  }
328  } else {
329  switch (type) {
330  case 'D':
331  value = &reader->singleDouble[index];
332  break;
333  case 'F':
334  value = &reader->singleFloat[index];
335  break;
336  case 'I':
337  value = &reader->singleInt[index];
338  break;
339  case 'B':
340  value = &reader->singleBool[index];
341  break;
342  }
343  }
344 
345  reader->tree->SetBranchAddress(name, value);
346  }
347 
349  if (multiple) {
350  switch (type) {
351  case 'D': {
352  const std::vector<Double_t> *values = static_cast<const std::vector<Double_t> *>(ptr);
353  if (!values)
354  values = &reader->multiDouble[index].second;
355  for (std::vector<Double_t>::const_iterator iter = values->begin(); iter != values->end(); iter++)
356  reader->values.add(name, *iter);
357  break;
358  }
359  case 'F': {
360  const std::vector<Float_t> *values = static_cast<const std::vector<Float_t> *>(ptr);
361  if (!values)
362  values = &reader->multiFloat[index].second;
363  for (std::vector<Float_t>::const_iterator iter = values->begin(); iter != values->end(); iter++)
364  reader->values.add(name, *iter);
365  break;
366  }
367  case 'I': {
368  const std::vector<Int_t> *values = static_cast<const std::vector<Int_t> *>(ptr);
369  if (!values)
370  values = &reader->multiInt[index].second;
371  for (std::vector<Int_t>::const_iterator iter = values->begin(); iter != values->end(); iter++)
372  reader->values.add(name, *iter);
373  break;
374  }
375  case 'B': {
376  const std::vector<Bool_t> *values = static_cast<const std::vector<Bool_t> *>(ptr);
377  if (!values)
378  values = &reader->multiBool[index].second;
379  for (std::vector<Bool_t>::const_iterator iter = values->begin(); iter != values->end(); iter++)
380  reader->values.add(name, *iter);
381  break;
382  }
383  }
384  } else {
385  double value = 0.0;
386 
387  switch (type) {
388  case 'D':
389  value = ptr ? *(const Double_t *)ptr : reader->singleDouble[index];
390  break;
391  case 'F':
392  value = ptr ? *(const Float_t *)ptr : reader->singleFloat[index];
393  break;
394  case 'I':
395  value = ptr ? *(const Int_t *)ptr : reader->singleInt[index];
396  break;
397  case 'B':
398  value = ptr ? *(const Bool_t *)ptr : reader->singleBool[index];
399  break;
400  }
401 
402  if (!optional || value != optVal)
403  reader->values.add(name, value);
404  }
405  }
406 
407 #define TREEREADER_ADD_IMPL(T, C) \
408  template <> \
409  void TreeReader::addSingle<T>(AtomicId name, const T *value, bool opt) { \
410  addTypeSingle(name, value, C, opt); \
411  } \
412  \
413  template <> \
414  void TreeReader::addMulti(AtomicId name, const std::vector<T> *value) { \
415  addTypeMulti(name, value, C); \
416  }
417 
418  TREEREADER_ADD_IMPL(Double_t, 'D')
419  TREEREADER_ADD_IMPL(Float_t, 'F')
422 
423 #undef TREEREADER_ADD_IMPL
424 
425 } // namespace PhysicsTools
type
Definition: HCALResponse.h:21
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
static std::pair< void *, std::vector< T > > makeMulti()
Definition: TreeReader.cc:169
Variable::ValueList fill()
Definition: TreeReader.cc:286
#define nullptr
std::vector< AtomicId > variables() const
Definition: TreeReader.cc:296
int Bool
Definition: Types.h:100
Variable::ValueList values
Definition: TreeReader.h:107
reader
Definition: DQM.py:105
void setTree(TTree *tree)
Definition: TreeReader.cc:55
#define TREEREADER_ADD_IMPL(T, C)
Definition: TreeReader.cc:407
Cheap generic unique keyword identifier class.
Definition: AtomicId.h:31
std::vector< std::pair< void *, std::vector< Float_t > > > multiFloat
Definition: TreeReader.h:97
std::vector< Bool > singleBool
Definition: TreeReader.h:104
Main interface class to the generic discriminator computer framework.
Definition: MVAComputer.h:39
std::vector< std::pair< void *, std::vector< Int_t > > > multiInt
Definition: TreeReader.h:98
std::vector< std::pair< void *, std::vector< Double_t > > > multiDouble
Definition: TreeReader.h:96
void automaticAdd(bool skipTarget=false, bool skipWeight=false)
Definition: TreeReader.cc:213
const std::complex< double > I
Definition: I.h:8
optional
Definition: Types.py:167
Definition: value.py:1
static const std::string B
Helper class that can contain an list of identifier-value pairs.
Definition: Variable.h:77
void setOptional(AtomicId name, bool opt, double optVal=kOptVal)
Definition: TreeReader.cc:117
void addTypeSingle(AtomicId name, const void *value, char type, bool opt)
Definition: TreeReader.cc:128
void update(TreeReader *reader) const
Definition: TreeReader.cc:304
uint64_t loop(const MVAComputer *mva)
Definition: TreeReader.cc:260
static const double kOptVal
Definition: TreeReader.h:54
unsigned long long uint64_t
Definition: Time.h:13
void addBranch(const std::string &expression, AtomicId name=AtomicId(), bool opt=true)
Definition: TreeReader.cc:60
std::map< AtomicId, Value > valueMap
Definition: TreeReader.h:106
std::vector< std::pair< void *, std::vector< Bool_t > > > multiBool
Definition: TreeReader.h:99
std::vector< Double_t > singleDouble
Definition: TreeReader.h:101
void fill(AtomicId name, TreeReader *reader) const
Definition: TreeReader.cc:348
std::vector< Float_t > singleFloat
Definition: TreeReader.h:102
Definition: tree.py:1
TreeReader & operator=(const TreeReader &orig)
Definition: TreeReader.cc:35
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
void addTypeMulti(AtomicId name, const void *value, char type)
Definition: TreeReader.cc:173
std::vector< Int_t > singleInt
Definition: TreeReader.h:103
void add(AtomicId id, double value)
Definition: Variable.h:103