CMS 3D CMS Logo

TableOutputBranches.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_NanoAOD_TableOutputBranches_h
2 #define PhysicsTools_NanoAOD_TableOutputBranches_h
3 
4 #include <string>
5 #include <vector>
6 #include <TTree.h>
11 
13 public:
16  if (desc->className() != "nanoaod::FlatTable")
17  throw cms::Exception("Configuration", "NanoAODOutputModule can only write out nanoaod::FlatTable objects");
18  }
19 
21  void branch(TTree &tree);
22 
25  void fill(const edm::OccurrenceForOutput &iWhatever, TTree &tree, bool extensions);
26 
27 private:
30  bool m_singleton = false;
33  UInt_t m_counter;
34  struct NamedBranchPtr {
36  TBranch *branch;
38  const std::string &atitle,
39  const std::string &rootType,
40  TBranch *branchptr = nullptr)
41  : name(aname), title(atitle), rootTypeCode(rootType), branch(branchptr) {}
42  };
43  TBranch *m_counterBranch = nullptr;
44  std::vector<NamedBranchPtr> m_floatBranches;
45  std::vector<NamedBranchPtr> m_intBranches;
46  std::vector<NamedBranchPtr> m_int8Branches;
47  std::vector<NamedBranchPtr> m_uint8Branches;
48  std::vector<NamedBranchPtr> m_uint32Branches;
49  std::vector<NamedBranchPtr> m_doubleBranches;
51 
52  template <typename T>
53  void fillColumn(NamedBranchPtr &pair, const nanoaod::FlatTable &tab) {
54  int idx = tab.columnIndex(pair.name);
55  if (idx == -1)
56  throw cms::Exception("LogicError", "Missing column in input for " + m_baseName + "_" + pair.name);
57  pair.branch->SetAddress(
58  tab.size() == 0 ? static_cast<T *>(nullptr)
59  : const_cast<T *>(&tab.columnData<T>(idx).front())); // SetAddress should take a const * !
60  }
61 };
62 
63 #endif
std::vector< NamedBranchPtr > m_uint32Branches
void defineBranchesFromFirstEvent(const nanoaod::FlatTable &tab)
std::vector< NamedBranchPtr > m_intBranches
std::vector< NamedBranchPtr > m_int8Branches
int columnIndex(const std::string &name) const
Definition: FlatTable.cc:3
auto columnData(unsigned int column) const
get a column by index (const)
Definition: FlatTable.h:73
enum TableOutputBranches::@871 m_extension
std::vector< NamedBranchPtr > m_uint8Branches
TableOutputBranches(const edm::BranchDescription *desc, const edm::EDGetToken &token)
void fill(const edm::OccurrenceForOutput &iWhatever, TTree &tree, bool extensions)
void branch(TTree &tree)
edm::EDGetToken m_token
NamedBranchPtr(const std::string &aname, const std::string &atitle, const std::string &rootType, TBranch *branchptr=nullptr)
Definition: tree.py:1
long double T
std::vector< NamedBranchPtr > m_floatBranches
unsigned int size() const
Definition: FlatTable.h:57
void fillColumn(NamedBranchPtr &pair, const nanoaod::FlatTable &tab)
std::vector< NamedBranchPtr > m_doubleBranches