ATLAS Offline Software
Loading...
Searching...
No Matches
HistogramFillerTree.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef AthenaMonitoringKernel_HistogramFiller_HistogramFillerTree_h
6#define AthenaMonitoringKernel_HistogramFiller_HistogramFillerTree_h
7
8#include "TTree.h"
9
10#include "HistogramFiller.h"
11#include <boost/algorithm/string.hpp>
12
13namespace Monitored {
14 template <typename T> void scalarFillerFunc(TBranch* branch, const IMonitoredVariable& var);
15 template <> void scalarFillerFunc<std::string>(TBranch* branch, const IMonitoredVariable& var);
16 template <typename T> void vectorFillerFunc(TBranch* branch, const IMonitoredVariable& var);
17 template <> void vectorFillerFunc<std::string>(TBranch* branch, const IMonitoredVariable& var);
18
23 public:
24 HistogramFillerTree(const HistogramDef& definition, std::shared_ptr<IHistogramProvider> provider)
25 : HistogramFiller(definition, std::move(provider)) {
27 }
28
29 virtual unsigned fill( const HistogramFiller::VariablesPack& vars ) const override {
30
31 if (vars.cut) {
32 const size_t maskSize = vars.cut->size();
33 // Abort if no cut entries or first (and only) entry is false
34 if (maskSize == 0 || (maskSize == 1 && !vars.cut->get(0))) { return 0; }
35
36 if (ATH_UNLIKELY(maskSize > 1)) {
37 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
38 log << MSG::WARNING << "HistogramFillerTree (" << m_histDef->alias
39 << ") does not support more than a single entry being filled at a time\n"
40 << "so a cut mask with > 1 entry doesn't make sense. Using first entry only." << endmsg;
41 if (!vars.cut->get(0)) { return 0; }
42 }
43 }
44
45 if (ATH_UNLIKELY(vars.size() != m_branchDefs.size())) {
46 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
47 log << MSG::ERROR << "Mismatch of passed variables and expected variables for " << m_histDef->alias
48 << "(" << vars.size() << ", " << m_branchDefs.size() << ")" << endmsg;
49 return 0;
50 }
51
52 auto tree = this->histogram<TTree>();
53 if (tree->GetListOfBranches()->GetEntries() == 0) {
55 }
56 auto branchList = tree->GetListOfBranches();
57 // following logic allows us to skip branches that were badly-defined
58 size_t idx = 0, idxgood = 0;
59 for (const auto& brdef : m_branchDefs) {
60 if (ATH_UNLIKELY(brdef.second == "IGNORE")) {
61 ++idx; continue;
62 }
63 TBranch* branch = static_cast<TBranch*>(branchList->At(idxgood));
64 m_fillerFunctions[idx](branch, *vars[idx]);
65 ++idx; ++idxgood;
66 }
67 for (Int_t i = 0; i < branchList->GetEntries(); ++i) {
68
69
70 }
71 tree->SetEntries(tree->GetEntries() + 1);
72 return 1;
73 }
74
75 private:
76 std::vector<std::pair<std::string, std::string>> m_branchDefs;
77 std::vector<std::function<void(TBranch*, const IMonitoredVariable&)>> m_fillerFunctions;
78
80 std::vector<std::string> tokenized;
81 boost::split(tokenized, m_histDef->treeDef, [](char c){ return c ==':'; });
82 for (const auto& token : tokenized) {
83 auto ipart = token.find('/');
84 if (ipart == std::string::npos) {
85 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
86 log << MSG::ERROR << "Tree " << m_histDef->alias << ": Badly formed variable definition " << token
87 << "; skipping" << endmsg;
88 continue;
89 }
90 auto branch = token.substr(0, ipart);
91 auto type = token.substr(ipart+1);
92 std::function<void(TBranch*, const IMonitoredVariable&)> fillerFunc;
93 // scalar
94 if (type == "B") fillerFunc = scalarFillerFunc<Char_t>;
95 else if (type == "b") fillerFunc = scalarFillerFunc<UChar_t>;
96 else if (type == "S") fillerFunc = scalarFillerFunc<Short_t>;
97 else if (type == "s") fillerFunc = scalarFillerFunc<UShort_t>;
98 else if (type == "I") fillerFunc = scalarFillerFunc<Int_t>;
99 else if (type == "i") fillerFunc = scalarFillerFunc<UInt_t>;
100 else if (type == "F") fillerFunc = scalarFillerFunc<Float_t>;
101 else if (type == "f") fillerFunc = scalarFillerFunc<Float16_t>;
102 else if (type == "D") fillerFunc = scalarFillerFunc<Double_t>;
103 else if (type == "d") fillerFunc = scalarFillerFunc<Double32_t>;
104 else if (type == "L") fillerFunc = scalarFillerFunc<Long64_t>;
105 else if (type == "l") fillerFunc = scalarFillerFunc<ULong64_t>;
106 else if (type == "O") fillerFunc = scalarFillerFunc<Bool_t>;
107 else if (type == "string") fillerFunc = scalarFillerFunc<std::string>;
108 else if (type == "vector<char>") fillerFunc = vectorFillerFunc<char>;
109 else if (type == "vector<unsigned char>") fillerFunc = vectorFillerFunc<unsigned char>;
110 else if (type == "vector<int>") fillerFunc = vectorFillerFunc<int>;
111 else if (type == "vector<unsigned int>") fillerFunc = vectorFillerFunc<unsigned int>;
112 else if (type == "vector<float>") fillerFunc = vectorFillerFunc<float>;
113 else if (type == "vector<double>") fillerFunc = vectorFillerFunc<double>;
114 else if (type == "vector<long>") fillerFunc = vectorFillerFunc<long>;
115 else if (type == "vector<unsigned long>") fillerFunc = vectorFillerFunc<unsigned long>;
116 else if (type == "vector<string>") fillerFunc = vectorFillerFunc<std::string>;
117 else if (type == "C") {
118 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
119 log << MSG::ERROR << "Tree " << m_histDef->alias << ": Branch type \"C\" not supported for branch" << branch
120 << "; please use \"string\"" << endmsg;
121 type = "IGNORE";
122 } else {
123 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
124 log << MSG::ERROR << "Tree " << m_histDef->alias << ": Unrecognized branch type " << type << " for branch " << branch
125 << "; ignoring branch " << endmsg;
126 type = "IGNORE";
127 }
128 m_branchDefs.emplace_back(branch, type);
129 m_fillerFunctions.push_back(std::move(fillerFunc));
130 }
131 }
132
133 template <typename T>
134 void branchHelper(TTree* tree, const std::pair<std::string, std::string>& brdef) const {
135 std::vector<T> dummy;
136 std::vector<T>* dummyptr = &dummy;
137 // we only need dummy object to exist through the end of this method
138 tree->Branch(brdef.first.c_str(), &dummyptr);
139 }
140
141 void createBranches(TTree* tree) const {
142 for (const auto& brdef : m_branchDefs) {
143 if (brdef.second == "vector<char>") branchHelper<char>(tree, brdef);
144 else if (brdef.second == "vector<unsigned char>") branchHelper<unsigned char>(tree, brdef);
145 else if (brdef.second == "vector<int>") branchHelper<int>(tree, brdef);
146 else if (brdef.second == "vector<unsigned int>") branchHelper<unsigned int>(tree, brdef);
147 else if (brdef.second == "vector<float>") branchHelper<float>(tree, brdef);
148 else if (brdef.second == "vector<double>") branchHelper<double>(tree, brdef);
149 else if (brdef.second == "vector<long>") branchHelper<long>(tree, brdef);
150 else if (brdef.second == "vector<unsigned long>") branchHelper<unsigned long>(tree, brdef);
151 else if (brdef.second == "vector<string>") branchHelper<std::string>(tree, brdef);
152 else if (brdef.second == "string") {
153 std::string dummy; tree->Branch(brdef.first.c_str(), &dummy);
154 } else if (brdef.second != "IGNORE") {
155 tree->Branch(brdef.first.c_str(), (void*) nullptr, (brdef.first+"/"+brdef.second).c_str());
156 }
157 }
158 }
159 };
160
161 // helper functions for filling branches
162 template <typename T>
163 void scalarFillerFunc(TBranch* branch, const IMonitoredVariable& var) {
164 if (ATH_UNLIKELY(var.size() == 0)) {
165 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
166 log << MSG::WARNING << "Tree " << branch->GetTree()->GetName() << ": Empty value passed to scalar branch fill for"
167 << branch->GetName() << endmsg;
168 return;
169 }
170 T tofill = var.get(0);
171 T* ptr = &tofill;
172 branch->SetAddress(ptr);
173 branch->Fill();
174 }
175
176 // specialization for string
177 template <>
178 void scalarFillerFunc<std::string>(TBranch* branch, const IMonitoredVariable& var) {
179 if (ATH_UNLIKELY(var.size() == 0)) {
180 MsgStream log(Athena::getMessageSvc(), "HistogramFillerTree");
181 log << MSG::WARNING << "Tree " << branch->GetTree()->GetName() << ": Empty value passed to scalar branch fill for"
182 << branch->GetName() << endmsg;
183 return;
184 }
185 std::string tofill{var.getString(0)};
186 branch->SetObject(&tofill);
187 //*static_cast<std::string*>(branch->GetObject()) = tofill;
188 branch->Fill();
189 }
190
191 template <typename T>
192 void vectorFillerFunc(TBranch* branch, const IMonitoredVariable& var) {
193 std::vector<T> tofill;
194 tofill.reserve(var.size());
195 for (size_t i = 0; i < var.size(); i++) {
196 tofill.push_back(var.get(i));
197 }
198 std::vector<T>* tofillptr = &tofill;
199 branch->SetAddress(&tofillptr);
200 branch->Fill();
201 }
202
203 // specialization for string
204 template <>
205 void vectorFillerFunc<std::string>(TBranch* branch, const IMonitoredVariable& var) {
206 std::vector<std::string> tofill;
207 tofill.reserve(var.size());
208 for (size_t i = 0; i < var.size(); i++) {
209 tofill.push_back(var.getString(i));
210 }
211 auto* tofillptr = &tofill;
212 branch->SetAddress(&tofillptr);
213 branch->Fill();
214 }
215}
216
217#endif // AthenaMonitoringKernel_HistogramFiller_HistogramFillerTree_h
#define endmsg
#define ATH_UNLIKELY(x)
std::vector< std::function< void(TBranch *, const IMonitoredVariable &)> > m_fillerFunctions
std::vector< std::pair< std::string, std::string > > m_branchDefs
virtual unsigned fill(const HistogramFiller::VariablesPack &vars) const override
Method that actually fills the ROOT object.
HistogramFillerTree(const HistogramDef &definition, std::shared_ptr< IHistogramProvider > provider)
void branchHelper(TTree *tree, const std::pair< std::string, std::string > &brdef) const
void createBranches(TTree *tree) const
HistogramFiller(const HistogramDef &histDef, std::shared_ptr< IHistogramProvider > histogramProvider)
Default constructor.
std::shared_ptr< HistogramDef > m_histDef
virtual double get(size_t) const =0
virtual size_t size() const =0
gives size of vector representation
IMessageSvc * getMessageSvc(bool quiet=false)
Generic monitoring tool for athena components.
void vectorFillerFunc< std::string >(TBranch *branch, const IMonitoredVariable &var)
void scalarFillerFunc< std::string >(TBranch *branch, const IMonitoredVariable &var)
void vectorFillerFunc(TBranch *branch, const IMonitoredVariable &var)
void scalarFillerFunc(TBranch *branch, const IMonitoredVariable &var)
STL namespace.
the internal class used to keep parsed Filler properties
helper class to pass variables to fillers
const Monitored::IMonitoredVariable * cut
pointer to cut mask variable, typically absent
size_t size() const
number of variables in the pack ( not counting the weight and mask )
TChain * tree