hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
HistoManager.cxx
Go to the documentation of this file.
1#include "HistoManager.h"
2#include <iostream>
3#include "TKey.h"
4#include "TClass.h"
5#include <fstream>
6#include <iomanip>
7#include <vector>
8
10 HistoManager("default");
11 m_name = "default";
12}
13
14HistoManager::HistoManager(const std::string& inputName) {
15 m_name = inputName;
16}
17
19
20 for (it1d it = histos1d.begin(); it!=histos1d.end(); ++it) {
21 if (it->second) {
22 delete (it->second);
23 (it->second) = nullptr;
24 }
25 }
26
27 histos1d.clear();
28
29 for (it2d it = histos2d.begin(); it!=histos2d.end(); ++it) {
30 if (it->second) {
31 delete (it->second);
32 (it->second) = nullptr;
33 }
34 }
35
36 histos2d.clear();
37
38 for (it3d it = histos3d.begin(); it!=histos3d.end(); ++it) {
39 if (it->second) {
40 delete (it->second);
41 (it->second) = nullptr;
42 }
43 }
44
45 histos3d.clear();
46
47}
48
50
52
53 if (debug_ > 0) std::cout << "[HistoManager] DefineHistos" << std::endl;
54 std::string h_name = "";
55 for (auto hist : _h_configs.items()) {
56 h_name = m_name+"_" + hist.key() ;
57 std::size_t found = (hist.key()).find_last_of("_");
58 std::string extension = hist.key().substr(found+1);
59
60 if (extension == "h") {
61 histos1d[h_name] = plot1D(h_name,hist.value().at("xtitle"),
62 hist.value().at("bins"),
63 hist.value().at("minX"),
64 hist.value().at("maxX"));
65
66 std::string ytitle = hist.value().at("ytitle");
67
68 histos1d[h_name]->GetYaxis()->SetTitle(ytitle.c_str());
69
70 if (hist.value().contains("labels")) {
71 std::vector<std::string> labels = hist.value().at("labels").get<std::vector<std::string> >();
72
73 if (labels.size() < hist.value().at("bins")) {
74 std::cout<<"Cannot apply labels to histogram:"<<h_name<<std::endl;
75 }
76 else {
77 for (int i = 1; i<=hist.value().at("bins");++i)
78 histos1d[h_name]->GetXaxis()->SetBinLabel(i,labels[i-1].c_str());
79 }//bins
80 }//labels
81 }//1D histo
82
83 else if (extension == "hh") {
84 histos2d[h_name] = plot2D(h_name,
85 hist.value().at("xtitle"),hist.value().at("binsX"),hist.value().at("minX"),hist.value().at("maxX"),
86 hist.value().at("ytitle"),hist.value().at("binsY"),hist.value().at("minY"),hist.value().at("maxY"));
87 }
88
89 //3D histogram
90 else if (extension == "hhh") {
91 histos3d[h_name] = plot3D(h_name,
92 hist.value().at("xtitle"),
93 hist.value().at("binsX"),
94 hist.value().at("minX"),
95 hist.value().at("maxX"),
96 hist.value().at("ytitle"),
97 hist.value().at("binsY"),
98 hist.value().at("minY"),
99 hist.value().at("maxY"),
100 hist.value().at("ztitle"),
101 hist.value().at("binsZ"),
102 hist.value().at("minZ"),
103 hist.value().at("maxZ"));
104 }
105 }//loop on config
106}
107
108//This DefineHistos method is only used if you want to make multiple differently named copies of a particular json file histogram configuration.
109//The json histo config entry key value must contain a special string tag that indicates whether or not to make copies.
110//Must also provide a list of names that will be appended to the histo config key
111void HistoManager::DefineHistos(std::vector<std::string> histoCopyNames, std::string makeCopyJsonTag){
112 if (debug_ > 0) std::cout << "[HistoManager] DefineHistos" << std::endl;
113 std::string h_name = "";
114 if (debug_ > 0){
115 for (auto hist : _h_configs.items()){
116 std::cout<<hist.key()<<std::endl;
117 }
118 }
119 //std::cout<<_h_configs.items()<<std::endl;
120 for (auto hist : _h_configs.items()) {
121 bool singleCopy = true;
122 std::cout << "hist copy list size " << histoCopyNames.size() << std::endl;
123 for(int i = 0; i < histoCopyNames.size(); i++){
124 h_name = m_name+"_" + hist.key() ;
125 if (histoCopyNames.size() > 1 && std::string(hist.key()).find(makeCopyJsonTag) != std::string::npos){
126 h_name = m_name+"_"+ histoCopyNames.at(i) + "_" + hist.key() ;
127 singleCopy = false;
128 }
129 std::size_t found = (hist.key()).find_last_of("_");
130 std::string extension = hist.key().substr(found+1);
131 if(debug_){
132 std::cout << "DefineHisto: " << h_name << std::endl;
133 std::cout << extension << hist.value().at("xtitle") << std::endl;
134 }
135 if (extension == "h") {
136 histos1d[h_name] = plot1D(h_name,hist.value().at("xtitle"),
137 hist.value().at("bins"),
138 hist.value().at("minX"),
139 hist.value().at("maxX"));
140
141 std::string ytitle = hist.value().at("ytitle");
142
143 histos1d[h_name]->GetYaxis()->SetTitle(ytitle.c_str());
144
145 if (hist.value().contains("labels")) {
146 std::vector<std::string> labels = hist.value().at("labels").get<std::vector<std::string> >();
147
148 if (labels.size() < hist.value().at("bins")) {
149 std::cout<<"Cannot apply labels to histogram:"<<h_name<<std::endl;
150 }
151 else {
152 for (int i = 1; i<=hist.value().at("bins");++i)
153 histos1d[h_name]->GetXaxis()->SetBinLabel(i,labels[i-1].c_str());
154 }//bins
155 }//labels
156 }//1D histo
157
158 else if (extension == "hh") {
159 histos2d[h_name] = plot2D(h_name,
160 hist.value().at("xtitle"),hist.value().at("binsX"),hist.value().at("minX"),hist.value().at("maxX"),
161 hist.value().at("ytitle"),hist.value().at("binsY"),hist.value().at("minY"),hist.value().at("maxY"));
162 }
163
164 if(singleCopy)
165 break;
166
167 }
168 }//loop on config
169}
170
171void HistoManager::GetHistosFromFile(TFile* inFile, const std::string& name, const std::string& folder) {
172
173 //Todo: use name as regular expression.
174 //Todo: use folder to choose a folder.
175 TIter next(inFile->GetListOfKeys());
176 TKey *key;
177 while ((key = (TKey*)next())) {
178 std::string classType = key->GetClassName();
179 if (classType.find("TH1")!=std::string::npos)
180 histos1d[key->GetName()] = (TH1F*) key->ReadObj();
181 if (classType.find("TH2")!=std::string::npos)
182 histos2d[key->GetName()] = (TH2F*) key->ReadObj();
183 if (classType.find("TH3")!=std::string::npos)
184 histos3d[key->GetName()] = (TH3F*) key->ReadObj();
185 }
186}
187
188
189TH1F* HistoManager::plot1D(const std::string& name,const std::string& xtitle, int nbinsX, float xmin, float xmax) {
190 TH1F* h= new TH1F(name.c_str(),name.c_str(),nbinsX,xmin,xmax);
191 h->GetXaxis()->SetTitle(xtitle.c_str());
192 h->Sumw2();
193 return h;
194}
195
196TH1F* HistoManager::plot1D(const std::string& name,const std::string& xtitle, int nbinsX, double* axisX) {
197 TH1F* h= new TH1F(name.c_str(),name.c_str(),nbinsX,axisX);
198 h->GetXaxis()->SetTitle(xtitle.c_str());
199 h->Sumw2();
200 return h;
201}
202
203TH2F* HistoManager::plot2D(std::string name,
204 std::string xtitle, int nbinsX, float xmin, float xmax,
205 std::string ytitle, int nbinsY, float ymin, float ymax) {
206
207 TH2F* h = new TH2F(name.c_str(),name.c_str(),
208 nbinsX,xmin,xmax,
209 nbinsY,ymin,ymax);
210 h->GetXaxis()->SetTitle(xtitle.c_str());
211 h->GetYaxis()->SetTitle(ytitle.c_str());
212 h->Sumw2();
213 return h;
214}
215
216TH2F* HistoManager::plot2D(std::string name,
217 std::string xtitle, int nbinsX, double* axisX,
218 std::string ytitle, int nbinsY, double* axisY) {
219
220 TH2F * h = new TH2F(name.c_str(),name.c_str(),
221 nbinsX,axisX,
222 nbinsY,axisY);
223 h->GetXaxis()->SetTitle(xtitle.c_str());
224 h->GetYaxis()->SetTitle(ytitle.c_str());
225 h->Sumw2();
226 return h;
227}
228
229TH2F* HistoManager::plot2D(std::string name,
230 std::string xtitle, int nbinsX, const double* axisX,
231 std::string ytitle, int nbinsY, const double* axisY) {
232
233 TH2F * h = new TH2F(name.c_str(),name.c_str(),
234 nbinsX,axisX,
235 nbinsY,axisY);
236 h->GetXaxis()->SetTitle(xtitle.c_str());
237 h->GetYaxis()->SetTitle(ytitle.c_str());
238 h->Sumw2();
239 return h;
240}
241
242
243
244TH2F* HistoManager::plot2D(std::string name,
245 std::string xtitle, int nbinsX, double* axisX,
246 std::string ytitle, int nbinsY, float ymin, float ymax) {
247
248 TH2F * h = new TH2F(name.c_str(),name.c_str(),
249 nbinsX,axisX,
250 nbinsY,ymin,ymax);
251 h->GetXaxis()->SetTitle(xtitle.c_str());
252 h->GetYaxis()->SetTitle(ytitle.c_str());
253 h->Sumw2();
254 return h;
255}
256
257
258TH3F* HistoManager::plot3D(std::string name,
259 std::string xtitle, int nbinsX, double* axisX,
260 std::string ytitle, int nbinsY, double* axisY,
261 std::string ztitle, int nbinsZ, double* axisZ) {
262
263
264 TH3F* h = new TH3F(name.c_str(),name.c_str(),
265 nbinsX,axisX,
266 nbinsY,axisY,
267 nbinsZ,axisZ);
268
269
270 h->GetXaxis()->SetTitle(xtitle.c_str());
271 h->GetYaxis()->SetTitle(ytitle.c_str());
272 h->GetZaxis()->SetTitle(ztitle.c_str());
273 h->Sumw2();
274 return h;
275}
276
277
278
279
280TH3F* HistoManager::plot3D(std::string name,
281 std::string xtitle, int nbinsX, float xmin, float xmax,
282 std::string ytitle, int nbinsY, float ymin, float ymax,
283 std::string ztitle, int nbinsZ, float zmin, float zmax) {
284
285
286 TH3F* h = new TH3F(name.c_str(),name.c_str(),
287 nbinsX,xmin,xmax,
288 nbinsY,ymin,ymax,
289 nbinsZ,zmin,zmax);
290
291 h->GetXaxis()->SetTitle(xtitle.c_str());
292 h->GetYaxis()->SetTitle(ytitle.c_str());
293 h->GetZaxis()->SetTitle(ztitle.c_str());
294 h->Sumw2();
295 return h;
296}
297
298
299
301
302 for (it3d it = histos3d.begin(); it!=histos3d.end(); ++it) {
303 if (!it->second){
304 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
305 continue;
306 }
307 it->second->Sumw2();
308 }
309
310 for (it2d it = histos2d.begin(); it!=histos2d.end(); ++it) {
311 if (!(it->second)) {
312 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
313 continue;
314 }
315 it->second->Sumw2();
316 }
317
318 for (it1d it = histos1d.begin(); it!=histos1d.end(); ++it) {
319 if (!it->second){
320 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
321 continue;
322 }
323 it->second->Sumw2();
324 }
325
326}
327
328void HistoManager::Fill2DHisto(const std::string& histoName,float valuex, float valuey, float weight) {
329 if (histos2d[m_name+"_"+histoName])
330 histos2d[m_name+"_"+histoName]->Fill(valuex,valuey,weight);
331 else {
333 if (doPrintWarnings_) {
335 std::cout<<"ERROR::Fill2DHisto Histogram not found! "<<m_name+"_"+histoName<<std::endl;
336 else {
337 std::cout<<"Fill2DHisto::Printed max number of warnings " << maxWarnings_ << ". Stop"<<std::endl;
338 doPrintWarnings_ = false;
339 }
340 }
341 }
342}
343
344void HistoManager::Fill3DHisto(const std::string& histoName,float valuex, float valuey, float valuez, float weight) {
345 if (histos3d[m_name+"_"+histoName])
346 histos3d[m_name+"_"+histoName]->Fill(valuex,valuey,valuez, weight);
347 else {
349 if (doPrintWarnings_) {
351 std::cout<<"ERROR::Fill3DHisto Histogram not found! "<<m_name+"_"+histoName<<std::endl;
352 else {
353 std::cout<<"Fill3DHisto::Printed max number of warnings " << maxWarnings_ << ". Stop"<<std::endl;
354 doPrintWarnings_ = false;
355 }
356 }
357 }
358}
359
360
361
362
363void HistoManager::Fill1DHisto(const std::string& histoName,float value, float weight) {
364 if (histos1d[m_name+"_"+histoName])
365 histos1d[m_name+"_"+histoName]->Fill(value,weight);
366 else {
368 if (doPrintWarnings_) {
370 std::cout<<"ERROR::Fill1DHisto Histogram not found! "<<m_name+"_"+histoName<<std::endl;
371 else {
372 std::cout<<"Fill1DHisto::Printed max number of warnings " << maxWarnings_ << ". Stop"<<std::endl;
373 doPrintWarnings_ = false;
374 }
375 }
376 }
377}
378
379
380void HistoManager::loadHistoConfig(const std::string histoConfigFile) {
381
382 std::ifstream i_file(histoConfigFile);
383 i_file >> _h_configs;
384 if (debug_) {
385 for (auto& el : _h_configs.items())
386 std::cout << el.key() << " : " << el.value() << "\n";
387 }
388 i_file.close();
389
390}
391
392
393
394
395void HistoManager::saveHistos(TFile* outF,std::string folder) {
396 if (outF) outF->cd();
397 TDirectory* dir{nullptr};
398 std::cout<<folder.c_str()<<std::endl;
399 if (!folder.empty()) {
400 dir = outF->mkdir(folder.c_str());
401 dir->cd();
402 }
403 for (it3d it = histos3d.begin(); it!=histos3d.end(); ++it) {
404 if (!it->second){
405 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
406 continue;
407 }
408 it->second->Write();
409 }
410 for (it2d it = histos2d.begin(); it!=histos2d.end(); ++it) {
411 if (!(it->second)) {
412 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
413 continue;
414 }
415 it->second->Write();
416 }
417
418 for (it1d it = histos1d.begin(); it!=histos1d.end(); ++it) {
419 if (!it->second){
420 std::cout<<it->first<<" Null ptr in saving.."<<std::endl;
421 continue;
422 }
423 it->second->Write();
424 }
425
426 //dir->Write();
427 //if (dir) {delete dir; dir=0;}
428
429 Clear();
430
431}
432
433
virtual ~HistoManager()
void Fill2DHisto(const std::string &histoName, float valuex, float valuey, float weight=1.)
description
std::map< std::string, TH2F * > histos2d
description
virtual void DefineHistos()
Definition of histograms from json config.
bool doPrintWarnings_
description
virtual void GetHistosFromFile(TFile *inFile, const std::string &name, const std::string &folder="")
Get histograms from input file.
TH2F * plot2D(std::string name, std::string xtitle, int nbinsX, float xmin, float xmax, std::string ytitle, int nbinsY, float ymin, float ymax)
description
std::map< std::string, TH3F * >::iterator it3d
description
HistoManager()
default constructor
std::map< std::string, TH1F * >::iterator it1d
description
void Fill3DHisto(const std::string &histoName, float valuex, float valuey, float valuez, float weight=1.)
description
int maxWarnings_
description
std::map< std::string, TH3F * > histos3d
description
std::map< std::string, TH1F * > histos1d
description
std::map< std::string, TH2F * >::iterator it2d
description
TH3F * plot3D(std::string name, std::string xtitle, int nbinsX, float xmin, float xmax, std::string ytitle, int nbinsY, float ymin, float ymax, std::string ztitle, int nbinsZ, float zmin, float zmax)
description
virtual void Clear()
description
void Fill1DHisto(const std::string &histoName, float value, float weight=1.)
description
virtual void sumw2()
description
TH1F * plot1D(const std::string &name, const std::string &xtitle, int nbinsX, float xmin, float xmax)
description
json _h_configs
description
virtual void loadHistoConfig(const std::string histoConfigFile)
load histogram config
std::string m_name
description
int printWarnings_
description
virtual void saveHistos(TFile *outF=nullptr, std::string folder="")
save histograms
bool debug_
description