AliPhysics  1909eaa (1909eaa)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
THistManager.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-2016, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 #include <cfloat>
16 #include <cstring>
17 #include <iostream> // for unit tests
18 #include <sstream>
19 #include <string>
20 #include <exception>
21 #include <vector>
22 #include <TArrayD.h>
23 #include <TAxis.h>
24 #include <TH1.h>
25 #include <TH2.h>
26 #include <TH3.h>
27 #include <THnSparse.h>
28 #include <THashList.h>
29 #include <TMath.h>
30 #include <TObjArray.h>
31 #include <TObjString.h>
32 #include <TProfile.h>
33 #include <TString.h>
34 
35 #include "TBinning.h"
36 #include "THistManager.h"
37 
41 
43  TNamed(),
44  fHistos(NULL),
45  fIsOwner(true)
46 {
47 }
48 
49 THistManager::THistManager(const char *name):
50  TNamed(name, Form("Histogram container %s", name)),
51  fHistos(NULL),
52  fIsOwner(true)
53 {
54  fHistos = new THashList();
55  fHistos->SetName(Form("histos%s", name));
56  fHistos->SetOwner();
57 }
58 
60  if(fHistos && fIsOwner) delete fHistos;
61 }
62 
63 THashList* THistManager::CreateHistoGroup(const char *groupname) {
64  // At first step check whether the group already exists.
65  THashList *foundgroup = FindGroup(groupname);
66  if(foundgroup){
67  // already existing, don't create it again
68  return foundgroup;
69  }
70  TString parentname = basename(groupname);
71  THashList *parentgroup = FindGroup(parentname);
72  if(!parentgroup) parentgroup = CreateHistoGroup(parentname);
73  THashList *childgroup = new THashList();
74  childgroup->SetName(histname(groupname));
75  childgroup->SetOwner();
76  parentgroup->Add(childgroup);
77  return childgroup;
78 }
79 
80 TH1* THistManager::CreateTH1(const char *name, const char *title, int nbins, double xmin, double xmax, Option_t *opt){
81  TString dirname(basename(name)), hname(histname(name));
82  THashList *parent(FindGroup(dirname));
83  if(!parent) parent = CreateHistoGroup(dirname);
84  if(parent->FindObject(hname)){
85  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
86  return 0;
87  }
88  TH1* h = new TH1D(hname, title, nbins, xmin, xmax);
89  TString optionstring(opt);
90  optionstring.ToLower();
91  if(optionstring.Contains("s"))
92  h->Sumw2();
93  parent->Add(h);
94  return h;
95 }
96 
97 TH1* THistManager::CreateTH1(const char *name, const char *title, int nbins, const double *xbins, Option_t *opt){
98  TString dirname(basename(name)), hname(histname(name));
99  THashList *parent(FindGroup(dirname));
100  if(!parent) parent = CreateHistoGroup(dirname);
101  if(parent->FindObject(hname)){
102  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
103  return 0;
104  }
105  TH1* h = new TH1D(hname, title, nbins, xbins);
106  TString optionstring(opt);
107  optionstring.ToLower();
108  if(optionstring.Contains("s"))
109  h->Sumw2();
110  parent->Add(h);
111  return h;
112 }
113 
114 TH1* THistManager::CreateTH1(const char *name, const char *title, const TArrayD &xbins, Option_t *opt){
115  TString dirname(basename(name)), hname(histname(name));
116  THashList *parent(FindGroup(dirname));
117  if(!parent) parent = CreateHistoGroup(dirname);
118  if(parent->FindObject(hname)){
119  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
120  return 0;
121  }
122  TH1* h = new TH1D(hname, title, xbins.GetSize()-1, xbins.GetArray());
123  TString optionstring(opt);
124  optionstring.ToLower();
125  if(optionstring.Contains("s"))
126  h->Sumw2();
127  parent->Add(h);
128  return h;
129 }
130 
131 TH1* THistManager::CreateTH1(const char *name, const char *title, const TBinning &xbin, Option_t *opt){
132  TArrayD myxbins;
133  try{
134  xbin.CreateBinEdges(myxbins);
135  } catch(std::exception &e){
136  Fatal("THistManager::CreateTH1", "Exception raised: %s", e.what());
137  }
138  return CreateTH1(name, title, myxbins, opt);
139 }
140 
141 TH2* THistManager::CreateTH2(const char *name, const char *title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, Option_t *opt){
142  TString dirname(basename(name)), hname(histname(name));
143  THashList *parent(FindGroup(dirname));
144  if(!parent) parent = CreateHistoGroup(dirname);
145  if(parent->FindObject(hname)){
146  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
147  return 0;
148  }
149  TH2* h = new TH2D(hname, title, nbinsx, xmin, xmax, nbinsy, ymin, ymax);
150  TString optionstring(opt);
151  optionstring.ToLower();
152  if(optionstring.Contains("s"))
153  h->Sumw2();
154  parent->Add(h);
155  return h;
156 }
157 
158 TH2* THistManager::CreateTH2(const char *name, const char *title, int nbinsx, const double *xbins, int nbinsy, const double *ybins, Option_t *opt){
159  TString dirname(basename(name)), hname(histname(name));
160  THashList *parent(FindGroup(dirname));
161  if(!parent) parent = CreateHistoGroup(dirname);
162  if(parent->FindObject(hname)){
163  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
164  return 0;
165  }
166  TH2* h = new TH2D(hname, title, nbinsx, xbins, nbinsy, ybins);
167  TString optionstring(opt);
168  optionstring.ToLower();
169  if(optionstring.Contains("s"))
170  h->Sumw2();
171  parent->Add(h);
172  return h;
173 }
174 
175 TH2* THistManager::CreateTH2(const char *name, const char *title, const TArrayD &xbins, const TArrayD &ybins, Option_t *opt){
176  TString dirname(basename(name)), hname(histname(name));
177  THashList *parent(FindGroup(dirname));
178  if(!parent) parent = CreateHistoGroup(dirname);
179  if(parent->FindObject(hname)){
180  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
181  return 0;
182  }
183  TH2* h = new TH2D(hname, title, xbins.GetSize() - 1, xbins.GetArray(), ybins.GetSize() - 1, ybins.GetArray());
184  TString optionstring(opt);
185  optionstring.ToLower();
186  if(optionstring.Contains("s"))
187  h->Sumw2();
188  parent->Add(h);
189  return h;
190 }
191 
192 TH2* THistManager::CreateTH2(const char *name, const char *title, const TBinning &xbins, const TBinning &ybins, Option_t *opt){
193  TArrayD myxbins, myybins;
194  try{
195  xbins.CreateBinEdges(myxbins);
196  } catch (std::exception &e) {
197  Fatal("THistManager::CreateTH2 (x-dir)", "Exception raised: %s", e.what());
198  }
199  try{
200  ybins.CreateBinEdges(myybins);
201  } catch (std::exception &e) {
202  Fatal("THistManager::CreateTH2 (y-dir)", "Exception raised: %s", e.what());
203  }
204 
205  return CreateTH2(name, title, myxbins, myybins, opt);
206 }
207 
208 TH3* THistManager::CreateTH3(const char* name, const char* title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, int nbinsz, double zmin, double zmax, Option_t *opt) {
209  TString dirname(basename(name)), hname(histname(name));
210  THashList *parent(FindGroup(dirname));
211  if(!parent) parent = CreateHistoGroup(dirname);
212  if(parent->FindObject(hname)){
213  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
214  return 0;
215  }
216  TH3* h = new TH3D(hname, title, nbinsx, xmin, xmax, nbinsy, ymin, ymax, nbinsz, zmin, zmax);
217  TString optionstring(opt);
218  optionstring.ToLower();
219  if(optionstring.Contains("s"))
220  h->Sumw2();
221  parent->Add(h);
222  return h;
223 }
224 
225 TH3* THistManager::CreateTH3(const char* name, const char* title, int nbinsx, const double* xbins, int nbinsy, const double* ybins, int nbinsz, const double* zbins, Option_t *opt) {
226  TString dirname(basename(name)), hname(histname(name));
227  THashList *parent(FindGroup(dirname));
228  if(!parent) parent = CreateHistoGroup(dirname);
229  if(parent->FindObject(hname)){
230  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
231  return 0;
232  }
233  TH3* h = new TH3D(hname, title, nbinsx, xbins, nbinsy, ybins, nbinsz, zbins);
234  TString optionstring(opt);
235  optionstring.ToLower();
236  if(optionstring.Contains("s"))
237  h->Sumw2();
238  parent->Add(h);
239  return h;
240 }
241 
242 TH3* THistManager::CreateTH3(const char* name, const char* title, const TArrayD& xbins, const TArrayD& ybins, const TArrayD& zbins, Option_t *opt) {
243  TString dirname(basename(name)), hname(histname(name));
244  THashList *parent(FindGroup(dirname));
245  if(!parent) parent = CreateHistoGroup(dirname);
246  if(parent->FindObject(hname)){
247  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
248  return 0;
249  }
250  TH3* h = new TH3D(hname, title, xbins.GetSize()-1, xbins.GetArray(), ybins.GetSize()-1, ybins.GetArray(), zbins.GetSize()-1, zbins.GetArray());
251  TString optionstring(opt);
252  optionstring.ToLower();
253  if(optionstring.Contains("s"))
254  h->Sumw2();
255  parent->Add(h);
256  return h;
257 }
258 
259 TH3* THistManager::CreateTH3(const char *name, const char *title, const TBinning &xbins, const TBinning &ybins, const TBinning &zbins, Option_t *opt){
260  TArrayD myxbins, myybins, myzbins;
261  try{
262  xbins.CreateBinEdges(myxbins);
263  } catch (std::exception &e) {
264  Fatal("THistManager::CreateTH2 (x-dir)", "Exception raised: %s", e.what());
265  }
266  try{
267  ybins.CreateBinEdges(myybins);
268  } catch (std::exception &e) {
269  Fatal("THistManager::CreateTH2 (y-dir)", "Exception raised: %s", e.what());
270  }
271  try{
272  zbins.CreateBinEdges(myzbins);
273  } catch (std::exception &e) {
274  Fatal("THistManager::CreateTH2 (z-dir)", "Exception raised: %s", e.what());
275  }
276 
277  return CreateTH3(name, title, myxbins, myybins, myzbins);
278 }
279 
280 THnSparse* THistManager::CreateTHnSparse(const char *name, const char *title, int ndim, const int *nbins, const double *min, const double *max, Option_t *opt) {
281  TString dirname(basename(name)), hname(histname(name));
282  THashList *parent(FindGroup(dirname));
283  if(!parent) parent = CreateHistoGroup(dirname);
284  if(parent->FindObject(hname)){
285  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
286  return 0;
287  }
288  THnSparse* h = new THnSparseD(hname, title, ndim, nbins, min, max);
289  TString optionstring(opt);
290  optionstring.ToLower();
291  if(optionstring.Contains("s"))
292  h->Sumw2();
293  parent->Add(h);
294  return h;
295 }
296 
297 THnSparse* THistManager::CreateTHnSparse(const char *name, const char *title, int ndim, const TAxis **axes, Option_t *opt) {
298  TString dirname(basename(name)), hname(histname(name));
299  THashList *parent(FindGroup(dirname));
300  if(!parent) parent = CreateHistoGroup(dirname);
301  if(parent->FindObject(hname)){
302  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
303  return 0;
304  }
305  TArrayD xmin(ndim), xmax(ndim);
306  TArrayI nbins(ndim);
307  for(int idim = 0; idim < ndim; ++idim){
308  const TAxis &myaxis = *(axes[idim]);
309  nbins[idim] = myaxis.GetNbins();
310  xmin[idim] = myaxis.GetXmin();
311  xmax[idim] = myaxis.GetXmax();
312  }
313  THnSparseD *hsparse = new THnSparseD(hname, title, ndim, nbins.GetArray(), xmin.GetArray(), xmax.GetArray());
314  for(int id = 0; id < ndim; ++id)
315  *(hsparse->GetAxis(id)) = *(axes[id]);
316  TString optionstring(opt);
317  optionstring.ToLower();
318  if(optionstring.Contains("s"))
319  hsparse->Sumw2();
320  parent->Add(hsparse);
321  return hsparse;
322 }
323 
324 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, double xmin, double xmax, Option_t *opt) {
325  TString dirname(basename(name)), hname(histname(name));
326  THashList *parent(FindGroup(dirname));
327  if(!parent) parent = CreateHistoGroup(dirname);
328  if(parent->FindObject(hname))
329  Fatal("THistManager::CreateTProfile", "Object %s already exists in group %s", hname.Data(), dirname.Data());
330  TProfile *hist = new TProfile(hname, title, nbinsX, xmin, xmax, opt);
331  parent->Add(hist);
332 }
333 
334 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, const double* xbins, Option_t *opt) {
335  TString dirname(basename(name)), hname(histname(name));
336  THashList *parent(FindGroup(dirname));
337  if(!parent) parent = CreateHistoGroup(dirname);
338  if(parent->FindObject(hname))
339  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
340  TProfile *hist = new TProfile(hname, title, nbinsX, xbins, opt);
341  parent->Add(hist);
342 }
343 
344 void THistManager::CreateTProfile(const char* name, const char* title, const TArrayD& xbins, Option_t *opt){
345  TString dirname(basename(name)), hname(histname(name));
346  THashList *parent(FindGroup(dirname));
347  if(!parent) parent = CreateHistoGroup(dirname);
348  if(parent->FindObject(hname))
349  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
350  TProfile *hist = new TProfile(hname.Data(), title, xbins.GetSize()-1, xbins.GetArray(), opt);
351  parent->Add(hist);
352 }
353 
354 void THistManager::CreateTProfile(const char *name, const char *title, const TBinning &xbins, Option_t *opt){
355  TArrayD myxbins;
356  try{
357  xbins.CreateBinEdges(myxbins);
358  } catch (std::exception &e){
359  Fatal("THistManager::CreateProfile", "Exception raised: %s", e.what());
360  }
361  CreateTProfile(name, title, myxbins, opt);
362 }
363 
364 void THistManager::SetObject(TObject * const o, const char *group) {
365  THashList *parent(FindGroup(group));
366  if(!parent) CreateHistoGroup(group);
367  if(parent->FindObject(o->GetName())){
368  Fatal("THistManager::SetObject", "Parent %s does not exist", strcmp(group, "/") ? group : "");
369  return;
370  }
371  if(!(dynamic_cast<THnBase *>(o) || dynamic_cast<TH1 *>(o))){
372  Fatal("THistManager::SetObject", "Object %s is not of a histogram type",o->GetName());
373  return;
374  }
375  fHistos->Add(o);
376 }
377 
378 void THistManager::FillTH1(const char *name, double x, double weight, Option_t *opt) {
379  TString dirname(basename(name)), hname(histname(name));
380  THashList *parent(FindGroup(dirname));
381  if(!parent){
382  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
383  return;
384  }
385  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname));
386  if(!hist){
387  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
388  return;
389  }
390  TString optionstring(opt);
391  if(optionstring.Contains("w")){
392  // use bin width as weight
393  Int_t bin = hist->GetXaxis()->FindBin(x);
394  // check if not overflow or underflow bin
395  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
396  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
397  }
398  hist->Fill(x, weight);
399 }
400 
401 void THistManager::FillTH1(const char *name, const char *label, double weight, Option_t *opt) {
402  TString dirname(basename(name)), hname(histname(name));
403  THashList *parent(FindGroup(dirname));
404  if(!parent){
405  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
406  return;
407  }
408  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname));
409  if(!hist){
410  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
411  return;
412  }
413  TString optionstring(opt);
414  if(optionstring.Contains("w")){
415  // use bin width as weight
416  // get bin for label
417  Int_t bin = hist->GetXaxis()->FindBin(label);
418  // check if not overflow or underflow bin
419  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
420  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
421  }
422  hist->Fill(label, weight);
423 }
424 
425 void THistManager::FillTH2(const char *name, double x, double y, double weight, Option_t *opt) {
426  TString dirname(basename(name)), hname(histname(name));
427  THashList *parent(FindGroup(dirname));
428  if(!parent){
429  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
430  return;
431  }
432  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname));
433  if(!hist){
434  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
435  return;
436  }
437  TString optstring(opt);
438  Double_t myweight = optstring.Contains("w") ? 1. : weight;
439  if(optstring.Contains("wx")){
440  Int_t binx = hist->GetXaxis()->FindBin(x);
441  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
442  }
443  if(optstring.Contains("wy")){
444  Int_t biny = hist->GetYaxis()->FindBin(y);
445  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
446  }
447  hist->Fill(x, y, myweight);
448 }
449 
450 void THistManager::FillTH2(const char *name, double *point, double weight, Option_t *opt) {
451  TString dirname(basename(name)), hname(histname(name));
452  THashList *parent(FindGroup(dirname));
453  if(!parent){
454  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
455  return;
456  }
457  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname));
458  if(!hist){
459  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
460  return;
461  }
462  TString optstring(opt);
463  Double_t myweight = optstring.Contains("w") ? 1. : weight;
464  if(optstring.Contains("wx")){
465  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
466  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
467  }
468  if(optstring.Contains("wy")){
469  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
470  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
471  }
472  hist->Fill(point[0], point[1], weight);
473 }
474 
475 void THistManager::FillTH3(const char* name, double x, double y, double z, double weight, Option_t *opt) {
476  TString dirname(basename(name)), hname(histname(name));
477  THashList *parent(FindGroup(dirname));
478  if(!parent){
479  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
480  return;
481  }
482  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname));
483  if(!hist){
484  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
485  return;
486  }
487  TString optstring(opt);
488  Double_t myweight = optstring.Contains("w") ? 1. : weight;
489  if(optstring.Contains("wx")){
490  Int_t binx = hist->GetXaxis()->FindBin(x);
491  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
492  }
493  if(optstring.Contains("wy")){
494  Int_t biny = hist->GetYaxis()->FindBin(y);
495  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
496  }
497  if(optstring.Contains("wz")){
498  Int_t binz = hist->GetZaxis()->FindBin(z);
499  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
500  }
501  hist->Fill(x, y, z, weight);
502 }
503 
504 void THistManager::FillTH3(const char* name, const double* point, double weight, Option_t *opt) {
505  TString dirname(basename(name)), hname(histname(name));
506  THashList *parent(FindGroup(dirname));
507  if(!parent){
508  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
509  return;
510  }
511  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname));
512  if(!hist){
513  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
514  return;
515  }
516  TString optstring(opt);
517  Double_t myweight = optstring.Contains("w") ? 1. : weight;
518  if(optstring.Contains("wx")){
519  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
520  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
521  }
522  if(optstring.Contains("wy")){
523  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
524  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
525  }
526  if(optstring.Contains("wz")){
527  Int_t binz = hist->GetZaxis()->FindBin(point[2]);
528  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
529  }
530  hist->Fill(point[0], point[1], point[2], weight);
531 }
532 
533 void THistManager::FillTHnSparse(const char *name, const double *x, double weight, Option_t *opt) {
534  TString dirname(basename(name)), hname(histname(name));
535  THashList *parent(FindGroup(dirname));
536  if(!parent){
537  Fatal("THistManager::FillTHnSparse", "Parent group %s does not exist", dirname.Data());
538  return;
539  }
540  THnSparseD *hist = dynamic_cast<THnSparseD *>(parent->FindObject(hname));
541  if(!hist){
542  Fatal("THistManager::FillTHnSparse", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
543  return;
544  }
545  TString optstring(opt);
546  Double_t myweight = optstring.Contains("w") ? 1. : weight;
547  for(Int_t iaxis = 0; iaxis < hist->GetNdimensions(); iaxis++){
548  std::stringstream weighthandler;
549  weighthandler << "w" << iaxis;
550  if(optstring.Contains(weighthandler.str().c_str())){
551  Int_t bin = hist->GetAxis(iaxis)->FindBin(x[iaxis]);
552  if(bin != 0 && bin != hist->GetAxis(iaxis)->GetNbins()) myweight *= hist->GetAxis(iaxis)->GetBinWidth(bin);
553  }
554  }
555 
556  hist->Fill(x, weight);
557 }
558 
559 void THistManager::FillProfile(const char* name, double x, double y, double weight){
560  TString dirname(basename(name)), hname(histname(name));
561  THashList *parent(FindGroup(dirname));
562  if(!parent)
563  Fatal("THistManager::FillTProfile", "Parent group %s does not exist", dirname.Data());
564  TProfile *hist = dynamic_cast<TProfile *>(parent->FindObject(hname));
565  if(!hist)
566  Fatal("THistManager::FillTProfile", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
567  hist->Fill(x, y, weight);
568 }
569 
570 TObject *THistManager::FindObject(const char *name) const {
571  TString dirname(basename(name)), hname(histname(name));
572  THashList *parent(FindGroup(dirname));
573  if(!parent) return NULL;
574  return parent->FindObject(hname);
575 }
576 
578  TString dirname(basename(obj->GetName())), hname(histname(obj->GetName()));
579  THashList *parent(FindGroup(dirname));
580  if(!parent) return NULL;
581  return parent->FindObject(hname);
582 }
583 
584 THashList *THistManager::FindGroup(const char *dirname) const {
585  if(!strlen(dirname) || !strcmp(dirname, "/")) return fHistos;
586  // recursive find - avoids tokenizing filename
587  THashList *parentlist = FindGroup(basename(dirname));
588  if(parentlist) return static_cast<THashList *>(parentlist->FindObject(histname(dirname)));
589  return nullptr;
590 }
591 
593  int index = path.Last('/');
594  if(index < 0) return ""; // no directory structure
595  return TString(path(0, index));
596 }
597 
599  int index = path.Last('/');
600  if(index < 0) return path; // no directory structure
601  return TString(path(index+1, path.Length() - (index+1)));
602 }
603 
609 
611  fkArray(hmgr),
612  fCurrentPos(),
613  fNext(),
614  fDirection(dir)
615 {}
616 
618  fkArray(ref.fkArray),
619  fCurrentPos(ref.fCurrentPos),
620  fNext(ref.fNext),
621  fDirection(ref.fDirection)
622 {}
623 
625  if(this != &ref){
626  fkArray = ref.fkArray;
627  fCurrentPos = ref.fCurrentPos;
628  fNext = ref.fNext;
629  fDirection = ref.fDirection;
630  }
631  return *this;
632 }
633 
635  return fCurrentPos == other.fCurrentPos;
636 }
637 
639  if(fDirection == kTHMIforward)
640  fCurrentPos++;
641  else
642  fCurrentPos--;
643  return *this;
644 }
645 
647  iterator tmp(*this);
648  operator++();
649  return tmp;
650 }
651 
653  if(fDirection == kTHMIforward)
654  fCurrentPos--;
655  else
656  fCurrentPos++;
657  return *this;
658 };
659 
661  iterator tmp(*this);
662  operator--();
663  return tmp;
664 };
665 
667  if(fCurrentPos >=0 && fCurrentPos < fkArray->GetListOfHistograms()->GetEntries())
668  return fkArray->GetListOfHistograms()->At(fCurrentPos);
669  return NULL;
670 }
671 
672 
678 
679 namespace TestTHistManager {
680 
681  int THistManagerTestSuite::TestBuildSimpleHistograms(){
682  THistManager testmgr("testmgr");
683 
684  // Create histogram of each type
685  testmgr.CreateTH1("Test1D", "Test Histogram 1D", 1, 0., 1.);
686  testmgr.CreateTH2("Test2D", "Test Histogram 2D", 2, 0., 2., 10., 0., 10);
687  testmgr.CreateTH3("Test3D", "Test Histogram 3D", 3, 2, 6., 10., 0., 10., 50., 0., 50.);
688  int nbins[3] = {3, 3, 3}; double min[3] = {0., 0., 0}, max[3] = {6, 9, 12}; // dimensions for THnSparseTest
689  testmgr.CreateTHnSparse("TestNSparse", "Test Histogram NSparse", 3, nbins, min, max);
690  testmgr.CreateTProfile("TestProfile", "Test TProfile", 1, 0., 1);
691 
692  // Find histograms in the list (evaluate test)
693  // Tell user why test has failed
694  Bool_t found(true);
695  if(!dynamic_cast<TH1 *>(testmgr.GetListOfHistograms()->FindObject("Test1D"))){
696  std::cout << "Not found: Test1D" << std::endl;
697  found = false;
698  }
699  if(!dynamic_cast<TH2 *>(testmgr.GetListOfHistograms()->FindObject("Test2D"))){
700  std::cout << "Not found: Test2D" << std::endl;
701  found = false;
702  }
703  if(!dynamic_cast<TH3 *>(testmgr.GetListOfHistograms()->FindObject("Test3D"))){
704  std::cout << "Not found: Test3D" << std::endl;
705  found = false;
706  }
707  if(!dynamic_cast<THnSparse *>(testmgr.GetListOfHistograms()->FindObject("TestNSparse"))){
708  std::cout << "Not found: TestNSparse" << std::endl;
709  found = false;
710  }
711  if(!dynamic_cast<TProfile *>(testmgr.GetListOfHistograms()->FindObject("TestProfile"))){
712  std::cout << "Not found: TestProfile" << std::endl;
713  found = false;
714  }
715 
716  return found ? 0 : 1;
717  }
718 
719  int THistManagerTestSuite::TestBuildGroupedHistograms(){
720  THistManager testmgr("testmgr");
721 
722  // Creating 3 groups
723  testmgr.CreateTH1("Group1/Test1", "Test Histogram 1 in group 1", 1, 0., 1.);
724  testmgr.CreateTH1("Group1/Test2", "Test Histogram 2 in group 1", 1, 0., 1.);
725  testmgr.CreateTH2("Group2/Test1", "Test Histogram 1 in group 2", 1, 0., 1., 1, 0., 1.);
726  testmgr.CreateTH2("Group2/Test2", "Test Histogram 2 in group 2", 1, 0., 1., 1, 0., 1.);
727  testmgr.CreateTH3("Group3/Test1", "Test Histogram 1 in group 3", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
728  testmgr.CreateTH3("Group3/Test2", "Test Histogram 2 in group 3", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
729  // Creating histogram in a subgroup
730  testmgr.CreateTProfile("Group4/Subgroup1/Test1", "Test histogram for subgroup handling", 1, 0., 1);
731 
732  // Evalutate test
733  // Tell user why test has failed
734  Bool_t found(true);
735  THashList *currentdir(nullptr);
736  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group1")))){
737  std::cout << "Not found: Group1" << std::endl;
738  found = false;
739  } else {
740  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test1"))){
741  std::cout << "Not found in Group1: Test1" << std::endl;
742  found = false;
743  }
744  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test2"))){
745  std::cout << "Not found in Group1: Test2" << std::endl;
746  found = false;
747  }
748  }
749  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group2")))){
750  std::cout << "Not found: Group2" << std::endl;
751  found = false;
752  } else {
753  if(!dynamic_cast<TH2 *>(currentdir->FindObject("Test1"))){
754  std::cout << "Not found in Group2: Test1" << std::endl;
755  found = false;
756  }
757  if(!dynamic_cast<TH2 *>(currentdir->FindObject("Test2"))){
758  std::cout << "Not found in Group2: Test2" << std::endl;
759  found = false;
760  }
761  }
762  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group3")))){
763  std::cout << "Not found: Group3" << std::endl;
764  found = false;
765  } else {
766  if(!static_cast<TH3 *>(currentdir->FindObject("Test1"))){
767  std::cout << "Not found in Group3: Test1" << std::endl;
768  found = false;
769  }
770  if(!static_cast<TH3 *>(currentdir->FindObject("Test2"))){
771  std::cout << "Not found in Group3: Test2" << std::endl;
772  found = false;
773  }
774  }
775  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group4")))){
776  std::cout << "Not found: Group4" << std::endl;
777  found = false;
778  } else {
779  if(!(currentdir = dynamic_cast<THashList *>(currentdir->FindObject("Subgroup1")))){
780  std::cout << "Not found in Group4: Subgroup1" << std::endl;
781  found = false;
782  } else {
783  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test1"))){
784  std::cout << "Not found in Subgroup1: Test1" << std::endl;
785  found = false;
786  }
787  }
788  }
789  return found ? 0 : 1;
790  }
791 
792  int THistManagerTestSuite::TestFillSimpleHistograms(){
793  THistManager testmgr("testmgr");
794 
795  testmgr.CreateTH1("Test1", "Test fill 1D histogram", 1, 0., 1.);
796  testmgr.CreateTH2("Test2", "Test fill 2D histogram", 1, 0., 1., 1, 0., 1.);
797  testmgr.CreateTH3("Test3", "Test fill 3D histogram", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
798  int nbins[4] = {1,1,1,1}; double min[4] = {0.,0.,0.,0.}, max[4] = {1.,1.,1.,1.};
799  testmgr.CreateTHnSparse("TestN", "Test Fill THnSparse", 4, nbins, min, max);
800  testmgr.CreateTProfile("TestProfile", "Test fill Profile histogram", 1, 0., 1.);
801 
802  double point[4] = {0.5, 0.5, 0.5, 0.5};
803  for(int i = 0; i < 100; i++){
804  testmgr.FillTH1("Test1", 0.5);
805  testmgr.FillTH2("Test2", 0.5, 0.5);
806  testmgr.FillTH3("Test3", 0.5, 0.5, 0.5);
807  testmgr.FillProfile("TestProfile", 0.5, 1.);
808  testmgr.FillTHnSparse("TestN", point);
809  }
810 
811  // Evalutate test
812  // tell user why test has failed
813  bool success(true);
814 
815  TH1 *test1 = dynamic_cast<TH1 *>(testmgr.GetListOfHistograms()->FindObject("Test1"));
816  if(test1){
817  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
818  std::cout << "Test1: Mismatch in values, expected 100, found " << test1->GetBinContent(1) << std::endl;
819  success = false;
820  }
821  } else {
822  std::cout << "Not found: Test1" << std::endl;
823  success = false;
824  }
825 
826  TH2 *test2 = dynamic_cast<TH2 *>(testmgr.GetListOfHistograms()->FindObject("Test2"));
827  if(test2){
828  if(TMath::Abs(test2->GetBinContent(1, 1) - 100) > DBL_EPSILON){
829  std::cout << "Test2: Mismatch in values, expected 100, found " << test2->GetBinContent(1,1) << std::endl;
830  success = false;
831  }
832  } else {
833  std::cout << "Not found: Test2" << std::endl;
834  success = false;
835  }
836 
837  TH3 *test3 = dynamic_cast<TH3 *>(testmgr.GetListOfHistograms()->FindObject("Test3"));
838  if(test3){
839  if(TMath::Abs(test3->GetBinContent(1, 1, 1) - 100) > DBL_EPSILON){
840  std::cout << "Test3: Mismatch in values, expected 100, found " << test3->GetBinContent(1,1,1) << std::endl;
841  success = false;
842  }
843  } else {
844  std::cout << "Not found: Test3" << std::endl;
845  success = false;
846  }
847 
848  THnSparse *testN = dynamic_cast<THnSparse *>(testmgr.GetListOfHistograms()->FindObject("TestN"));
849  if(testN){
850  int index[4] = {1,1,1,1};
851  if(TMath::Abs(testN->GetBinContent(index) - 100) > DBL_EPSILON){
852  std::cout << "TestN: Mismatch in values, expected 100, found " << testN->GetBinContent(index) << std::endl;
853  success = false;
854  }
855  } else {
856  std::cout << "Not found: TestN" << std::endl;
857  success = false;
858  }
859 
860  TProfile *testProfile = dynamic_cast<TProfile *>(testmgr.GetListOfHistograms()->FindObject("TestProfile"));
861  if(testProfile){
862  if(TMath::Abs(testProfile->GetBinContent(1) - 1) > DBL_EPSILON){
863  std::cout << "TestProfile: Mismatch in values, expected 1, found " << testProfile->GetBinContent(1) << std::endl;
864  success = false;
865  }
866  } else {
867  std::cout << "Not found: TestProfile" << std::endl;
868  success = false;
869  }
870 
871  return success ? 0 : 1;
872  }
873 
874 
875  int THistManagerTestSuite::TestFillGroupedHistograms(){
876  THistManager testmgr("testmgr");
877 
878  // Creating 3 groups, 1 with 1D and 1 with 2D histograms, and a third with a subgroup and a TProfile
879  testmgr.CreateTH1("Group1/Test1", "Test 1 Group 1D", 1, 0., 1.);
880  testmgr.CreateTH1("Group1/Test2", "Test 2 Group 1D", 1, 0., 1.);
881  testmgr.CreateTH2("Group2/Test1", "Test 1 Group 2D", 1, 0., 1., 1, 0., 1.);
882  testmgr.CreateTH2("Group2/Test2", "Test 2 Group 2D", 1, 0., 1., 1, 0., 1.);
883  testmgr.CreateTProfile("Group3/Subgroup1/Test1", "Test 1 with subgroup", 1, 0., 1.);
884 
885  for(int i = 0; i < 100; i++){
886  testmgr.FillTH1("Group1/Test1", 0.5);
887  testmgr.FillTH1("Group1/Test2", 0.5);
888  testmgr.FillTH2("Group2/Test1", 0.5, 0.5);
889  testmgr.FillTH2("Group2/Test2", 0.5, 0.5);
890  testmgr.FillProfile("Group3/Subgroup1/Test1", 0.5, 1);
891  }
892 
893  // Evaluate test
894  bool success(true);
895 
896  THashList *currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group1"));
897  if(currentdir){
898  TH1 *test1 = dynamic_cast<TH1 *>(currentdir->FindObject("Test1"));
899  if(test1){
900  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
901  std::cout << "Group1/Test1: Value mismatch: expected 100, found " << test1->GetBinContent(1) << std::endl;
902  success = false;
903  }
904  } else {
905  std::cout << "Not found in Group1: Test1" << std::endl;
906  success = false;
907  }
908  test1 = dynamic_cast<TH1 *>(currentdir->FindObject("Test2"));
909  if(test1){
910  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
911  std::cout << "Group1/Test2: Value mismatch: expected 100, found " << test1->GetBinContent(1) << std::endl;
912  success = false;
913  }
914  } else {
915  std::cout << "Not found in Group1: Test2" << std::endl;
916  success = false;
917  }
918  } else {
919  std::cout << "Not found: Group1" << std::endl;
920  success = false;
921  }
922 
923  currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group2"));
924  if(currentdir){
925  TH2 *test2 = dynamic_cast<TH2 *>(currentdir->FindObject("Test1"));
926  if(test2){
927  if(TMath::Abs(test2->GetBinContent(1,1) - 100) > DBL_EPSILON){
928  std::cout << "Group2/Test1: Value mismatch: expected 100, found " << test2->GetBinContent(1,1) << std::endl;
929  success = false;
930  }
931  } else {
932  std::cout << "Not found in Group2: Test1" << std::endl;
933  success = false;
934  }
935  test2 = dynamic_cast<TH2 *>(currentdir->FindObject("Test2"));
936  if(test2){
937  if(TMath::Abs(test2->GetBinContent(1,1) - 100) > DBL_EPSILON){
938  std::cout << "Group2/Test2: Value mismatch: expected 100, found " << test2->GetBinContent(1,1) << std::endl;
939  success = false;
940  }
941  } else {
942  std::cout << "Not found in Group2: Test2" << std::endl;
943  success = false;
944  }
945  } else {
946  std::cout << "Not found: Group2" << std::endl;
947  success = false;
948  }
949 
950  currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group3"));
951  if(currentdir){
952  currentdir = dynamic_cast<THashList *>(currentdir->FindObject("Subgroup1"));
953  if(currentdir){
954  TProfile *testprofile = dynamic_cast<TProfile *>(currentdir->FindObject("Test1"));
955  if(testprofile){
956  if(TMath::Abs(testprofile->GetBinContent(1) - 1) > DBL_EPSILON){
957  std::cout << "Group3/Subgroup1/Test1: Value mismatch: expected 1, found " << testprofile->GetBinContent(1) << std::endl;
958  success = false;
959  }
960  } else {
961  std::cout << "Not found in Group3/Subgroup1: Test1" << std::endl;
962  success = false;
963  }
964  } else {
965  std::cout << "Not found in Group3: Subgroup1" << std::endl;
966  success = false;
967  }
968  } else {
969  std::cout << "Not found: Group3" << std::endl;
970  success = false;
971  }
972  return success ? 0 : 1;
973  }
974 
975  int TestRunAll(){
976  int testresult(0);
977  THistManagerTestSuite testsuite;
978 
979  std::cout << "Running test: Build simple" << std::endl;
980  testresult += testsuite.TestBuildSimpleHistograms();
981  std::cout << "Result after test: " << testresult << std::endl;
982 
983  std::cout << "Running test: Build grouped" << std::endl;
984  testresult += testsuite.TestBuildGroupedHistograms();
985  std::cout << "Result after test: " << testresult << std::endl;
986 
987  std::cout << "Running test: Fill Simple" << std::endl;
988  testresult += testsuite.TestFillSimpleHistograms();
989  std::cout << "Result after test: " << testresult << std::endl;
990 
991  std::cout << "Running test: Fill Grouped" << std::endl;
992  testresult += testsuite.TestFillGroupedHistograms();
993  std::cout << "Result after test: " << testresult << std::endl;
994 
995  return testresult;
996  }
997 
999  THistManagerTestSuite testsuite;
1000  return testsuite.TestBuildSimpleHistograms();
1001  }
1002 
1004  THistManagerTestSuite testsuite;
1005  return testsuite.TestBuildGroupedHistograms();
1006  }
1007 
1009  THistManagerTestSuite testsuite;
1010  return testsuite.TestFillSimpleHistograms();
1011  }
1012 
1014  THistManagerTestSuite testsuite;
1015  return testsuite.TestFillGroupedHistograms();
1016  }
1017 }
iterator & operator=(const iterator &rhs)
THashList * CreateHistoGroup(const char *groupname)
const Double_t ymax
Definition: AddTaskCFDStar.C:7
double Double_t
Definition: External.C:58
const char * title
Definition: MakeQAPdf.C:26
Definition: External.C:244
bool fIsOwner
Set the ownership.
Definition: THistManager.h:549
Int_t nbinsy
THashList * FindGroup(const char *dirname) const
void FillTH2(const char *hname, double x, double y, double weight=1., Option_t *opt="")
void FillTH3(const char *hname, double x, double y, double z, double weight=1., Option_t *opt="")
virtual void CreateBinEdges(TArrayD &binedges) const =0
const THistManager * fkArray
Definition: THistManager.h:156
void SetObject(TObject *const o, const char *group="/")
Interface for binnings used by the histogram handler.
Definition: TBinning.h:21
THMIDirection_t fDirection
Definition: THistManager.h:159
Bool_t operator!=(const iterator &aIter) const
void FillTHnSparse(const char *name, const double *x, double weight=1., Option_t *opt="")
TH2 * CreateTH2(const char *name, const char *title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, Option_t *opt="")
TObject * FindObject(const char *name) const
int Int_t
Definition: External.C:63
void CreateTProfile(const char *name, const char *title, int nbinsX, double xmin, double xmax, Option_t *opt="")
THashList * GetListOfHistograms() const
Definition: THistManager.h:504
Definition: External.C:252
TH1 * CreateTH1(const char *name, const char *title, int nbins, double xmin, double xmax, Option_t *opt="")
const Double_t zmin
Definition: External.C:228
Definition: External.C:212
void FillProfile(const char *name, double x, double y, double weight=1.)
THMIDirection_t
Direction for the iteration.
Definition: THistManager.h:72
stl-iterator for the histogram manager
Definition: THistManager.h:59
void FillTH1(const char *hname, double x, double weight=1., Option_t *opt="")
Int_t nbinsx
TString basename(const TString &path) const
Collection of tests for the THistManager.
Definition: THistManager.h:592
TString histname(const TString &path) const
Definition: External.C:220
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
TObject * operator*() const
const Double_t ymin
Definition: AddTaskCFDStar.C:6
Container class for histograms.
Definition: THistManager.h:43
const char Option_t
Definition: External.C:48
const Double_t zmax
const Int_t nbins
bool Bool_t
Definition: External.C:53
THnSparse * CreateTHnSparse(const char *name, const char *title, int ndim, const int *nbins, const double *min, const double *max, Option_t *opt="")
THashList * fHistos
List of histograms.
Definition: THistManager.h:548
Definition: External.C:196
TH3 * CreateTH3(const char *name, const char *title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, int nbinsz, double zmin, double zmax, Option_t *opt="")