AliPhysics  a5cd6b6 (a5cd6b6)
 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 THnSparse* THistManager::CreateTHnSparse(const char *name, const char *title, int ndim, const TBinning **axes, 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::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
330  return 0;
331  }
332  TArrayD xmin(ndim), xmax(ndim);
333  TArrayI nbins(ndim);
334  std::vector<TArrayD> binnings;
335  for(int idim = 0; idim < ndim; ++idim){
336  const TBinning &myaxis = *(axes[idim]);
337  TArrayD binEdges;
338  myaxis.CreateBinEdges(binEdges);
339  nbins[idim] = binEdges.GetSize()-1;
340  xmin[idim] = binEdges[0];
341  xmax[idim] = binEdges[binEdges.GetSize()-1];
342  binnings.push_back(binEdges);
343  }
344  THnSparseD *hsparse = new THnSparseD(hname, title, ndim, nbins.GetArray(), xmin.GetArray(), xmax.GetArray());
345  for(int id = 0; id < ndim; ++id){
346  TArrayD &binEdges = binnings[id];
347  hsparse->GetAxis(id)->Set(binEdges.GetSize()-1, binEdges.GetArray());
348  }
349  TString optionstring(opt);
350  optionstring.ToLower();
351  if(optionstring.Contains("s"))
352  hsparse->Sumw2();
353  parent->Add(hsparse);
354  return hsparse;
355 }
356 
357 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, double xmin, double xmax, Option_t *opt) {
358  TString dirname(basename(name)), hname(histname(name));
359  THashList *parent(FindGroup(dirname));
360  if(!parent) parent = CreateHistoGroup(dirname);
361  if(parent->FindObject(hname))
362  Fatal("THistManager::CreateTProfile", "Object %s already exists in group %s", hname.Data(), dirname.Data());
363  TProfile *hist = new TProfile(hname, title, nbinsX, xmin, xmax, opt);
364  parent->Add(hist);
365 }
366 
367 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, const double* xbins, Option_t *opt) {
368  TString dirname(basename(name)), hname(histname(name));
369  THashList *parent(FindGroup(dirname));
370  if(!parent) parent = CreateHistoGroup(dirname);
371  if(parent->FindObject(hname))
372  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
373  TProfile *hist = new TProfile(hname, title, nbinsX, xbins, opt);
374  parent->Add(hist);
375 }
376 
377 void THistManager::CreateTProfile(const char* name, const char* title, const TArrayD& xbins, Option_t *opt){
378  TString dirname(basename(name)), hname(histname(name));
379  THashList *parent(FindGroup(dirname));
380  if(!parent) parent = CreateHistoGroup(dirname);
381  if(parent->FindObject(hname))
382  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
383  TProfile *hist = new TProfile(hname.Data(), title, xbins.GetSize()-1, xbins.GetArray(), opt);
384  parent->Add(hist);
385 }
386 
387 void THistManager::CreateTProfile(const char *name, const char *title, const TBinning &xbins, Option_t *opt){
388  TArrayD myxbins;
389  try{
390  xbins.CreateBinEdges(myxbins);
391  } catch (std::exception &e){
392  Fatal("THistManager::CreateProfile", "Exception raised: %s", e.what());
393  }
394  CreateTProfile(name, title, myxbins, opt);
395 }
396 
397 void THistManager::SetObject(TObject * const o, const char *group) {
398  THashList *parent(FindGroup(group));
399  if(!parent) CreateHistoGroup(group);
400  if(parent->FindObject(o->GetName())){
401  Fatal("THistManager::SetObject", "Parent %s does not exist", strcmp(group, "/") ? group : "");
402  return;
403  }
404  if(!(dynamic_cast<THnBase *>(o) || dynamic_cast<TH1 *>(o))){
405  Fatal("THistManager::SetObject", "Object %s is not of a histogram type",o->GetName());
406  return;
407  }
408  fHistos->Add(o);
409 }
410 
411 void THistManager::FillTH1(const char *name, double x, double weight, Option_t *opt) {
412  TString dirname(basename(name)), hname(histname(name));
413  THashList *parent(FindGroup(dirname));
414  if(!parent){
415  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
416  return;
417  }
418  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname));
419  if(!hist){
420  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
421  return;
422  }
423  TString optionstring(opt);
424  if(optionstring.Contains("w")){
425  // use bin width as weight
426  Int_t bin = hist->GetXaxis()->FindBin(x);
427  // check if not overflow or underflow bin
428  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
429  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
430  }
431  hist->Fill(x, weight);
432 }
433 
434 void THistManager::FillTH1(const char *name, const char *label, double weight, Option_t *opt) {
435  TString dirname(basename(name)), hname(histname(name));
436  THashList *parent(FindGroup(dirname));
437  if(!parent){
438  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
439  return;
440  }
441  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname));
442  if(!hist){
443  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
444  return;
445  }
446  TString optionstring(opt);
447  if(optionstring.Contains("w")){
448  // use bin width as weight
449  // get bin for label
450  Int_t bin = hist->GetXaxis()->FindBin(label);
451  // check if not overflow or underflow bin
452  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
453  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
454  }
455  hist->Fill(label, weight);
456 }
457 
458 void THistManager::FillTH2(const char *name, double x, double y, double weight, Option_t *opt) {
459  TString dirname(basename(name)), hname(histname(name));
460  THashList *parent(FindGroup(dirname));
461  if(!parent){
462  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
463  return;
464  }
465  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname));
466  if(!hist){
467  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
468  return;
469  }
470  TString optstring(opt);
471  Double_t myweight = optstring.Contains("w") ? 1. : weight;
472  if(optstring.Contains("wx")){
473  Int_t binx = hist->GetXaxis()->FindBin(x);
474  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
475  }
476  if(optstring.Contains("wy")){
477  Int_t biny = hist->GetYaxis()->FindBin(y);
478  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
479  }
480  hist->Fill(x, y, myweight);
481 }
482 
483 void THistManager::FillTH2(const char *name, double *point, double weight, Option_t *opt) {
484  TString dirname(basename(name)), hname(histname(name));
485  THashList *parent(FindGroup(dirname));
486  if(!parent){
487  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
488  return;
489  }
490  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname));
491  if(!hist){
492  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
493  return;
494  }
495  TString optstring(opt);
496  Double_t myweight = optstring.Contains("w") ? 1. : weight;
497  if(optstring.Contains("wx")){
498  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
499  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
500  }
501  if(optstring.Contains("wy")){
502  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
503  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
504  }
505  hist->Fill(point[0], point[1], weight);
506 }
507 
508 void THistManager::FillTH3(const char* name, double x, double y, double z, double weight, Option_t *opt) {
509  TString dirname(basename(name)), hname(histname(name));
510  THashList *parent(FindGroup(dirname));
511  if(!parent){
512  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
513  return;
514  }
515  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname));
516  if(!hist){
517  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
518  return;
519  }
520  TString optstring(opt);
521  Double_t myweight = optstring.Contains("w") ? 1. : weight;
522  if(optstring.Contains("wx")){
523  Int_t binx = hist->GetXaxis()->FindBin(x);
524  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
525  }
526  if(optstring.Contains("wy")){
527  Int_t biny = hist->GetYaxis()->FindBin(y);
528  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
529  }
530  if(optstring.Contains("wz")){
531  Int_t binz = hist->GetZaxis()->FindBin(z);
532  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
533  }
534  hist->Fill(x, y, z, weight);
535 }
536 
537 void THistManager::FillTH3(const char* name, const double* point, double weight, Option_t *opt) {
538  TString dirname(basename(name)), hname(histname(name));
539  THashList *parent(FindGroup(dirname));
540  if(!parent){
541  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
542  return;
543  }
544  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname));
545  if(!hist){
546  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
547  return;
548  }
549  TString optstring(opt);
550  Double_t myweight = optstring.Contains("w") ? 1. : weight;
551  if(optstring.Contains("wx")){
552  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
553  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
554  }
555  if(optstring.Contains("wy")){
556  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
557  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
558  }
559  if(optstring.Contains("wz")){
560  Int_t binz = hist->GetZaxis()->FindBin(point[2]);
561  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
562  }
563  hist->Fill(point[0], point[1], point[2], weight);
564 }
565 
566 void THistManager::FillTHnSparse(const char *name, const double *x, double weight, Option_t *opt) {
567  TString dirname(basename(name)), hname(histname(name));
568  THashList *parent(FindGroup(dirname));
569  if(!parent){
570  Fatal("THistManager::FillTHnSparse", "Parent group %s does not exist", dirname.Data());
571  return;
572  }
573  THnSparseD *hist = dynamic_cast<THnSparseD *>(parent->FindObject(hname));
574  if(!hist){
575  Fatal("THistManager::FillTHnSparse", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
576  return;
577  }
578  TString optstring(opt);
579  Double_t myweight = optstring.Contains("w") ? 1. : weight;
580  for(Int_t iaxis = 0; iaxis < hist->GetNdimensions(); iaxis++){
581  std::stringstream weighthandler;
582  weighthandler << "w" << iaxis;
583  if(optstring.Contains(weighthandler.str().c_str())){
584  Int_t bin = hist->GetAxis(iaxis)->FindBin(x[iaxis]);
585  if(bin != 0 && bin != hist->GetAxis(iaxis)->GetNbins()) myweight *= hist->GetAxis(iaxis)->GetBinWidth(bin);
586  }
587  }
588 
589  hist->Fill(x, weight);
590 }
591 
592 void THistManager::FillProfile(const char* name, double x, double y, double weight){
593  TString dirname(basename(name)), hname(histname(name));
594  THashList *parent(FindGroup(dirname));
595  if(!parent)
596  Fatal("THistManager::FillTProfile", "Parent group %s does not exist", dirname.Data());
597  TProfile *hist = dynamic_cast<TProfile *>(parent->FindObject(hname));
598  if(!hist)
599  Fatal("THistManager::FillTProfile", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
600  hist->Fill(x, y, weight);
601 }
602 
603 TObject *THistManager::FindObject(const char *name) const {
604  TString dirname(basename(name)), hname(histname(name));
605  THashList *parent(FindGroup(dirname));
606  if(!parent) return NULL;
607  return parent->FindObject(hname);
608 }
609 
611  TString dirname(basename(obj->GetName())), hname(histname(obj->GetName()));
612  THashList *parent(FindGroup(dirname));
613  if(!parent) return NULL;
614  return parent->FindObject(hname);
615 }
616 
617 THashList *THistManager::FindGroup(const char *dirname) const {
618  if(!strlen(dirname) || !strcmp(dirname, "/")) return fHistos;
619  // recursive find - avoids tokenizing filename
620  THashList *parentlist = FindGroup(basename(dirname));
621  if(parentlist) return static_cast<THashList *>(parentlist->FindObject(histname(dirname)));
622  return nullptr;
623 }
624 
626  int index = path.Last('/');
627  if(index < 0) return ""; // no directory structure
628  return TString(path(0, index));
629 }
630 
632  int index = path.Last('/');
633  if(index < 0) return path; // no directory structure
634  return TString(path(index+1, path.Length() - (index+1)));
635 }
636 
642 
644  fkArray(hmgr),
645  fCurrentPos(),
646  fNext(),
647  fDirection(dir)
648 {}
649 
651  fkArray(ref.fkArray),
652  fCurrentPos(ref.fCurrentPos),
653  fNext(ref.fNext),
654  fDirection(ref.fDirection)
655 {}
656 
658  if(this != &ref){
659  fkArray = ref.fkArray;
660  fCurrentPos = ref.fCurrentPos;
661  fNext = ref.fNext;
662  fDirection = ref.fDirection;
663  }
664  return *this;
665 }
666 
668  return fCurrentPos == other.fCurrentPos;
669 }
670 
672  if(fDirection == kTHMIforward)
673  fCurrentPos++;
674  else
675  fCurrentPos--;
676  return *this;
677 }
678 
680  iterator tmp(*this);
681  operator++();
682  return tmp;
683 }
684 
686  if(fDirection == kTHMIforward)
687  fCurrentPos--;
688  else
689  fCurrentPos++;
690  return *this;
691 };
692 
694  iterator tmp(*this);
695  operator--();
696  return tmp;
697 };
698 
700  if(fCurrentPos >=0 && fCurrentPos < fkArray->GetListOfHistograms()->GetEntries())
701  return fkArray->GetListOfHistograms()->At(fCurrentPos);
702  return NULL;
703 }
704 
705 
711 
712 namespace TestTHistManager {
713 
714  int THistManagerTestSuite::TestBuildSimpleHistograms(){
715  THistManager testmgr("testmgr");
716 
717  // Create histogram of each type
718  testmgr.CreateTH1("Test1D", "Test Histogram 1D", 1, 0., 1.);
719  testmgr.CreateTH2("Test2D", "Test Histogram 2D", 2, 0., 2., 10., 0., 10);
720  testmgr.CreateTH3("Test3D", "Test Histogram 3D", 3, 2, 6., 10., 0., 10., 50., 0., 50.);
721  int nbins[3] = {3, 3, 3}; double min[3] = {0., 0., 0}, max[3] = {6, 9, 12}; // dimensions for THnSparseTest
722  testmgr.CreateTHnSparse("TestNSparse", "Test Histogram NSparse", 3, nbins, min, max);
723  testmgr.CreateTProfile("TestProfile", "Test TProfile", 1, 0., 1);
724 
725  // Find histograms in the list (evaluate test)
726  // Tell user why test has failed
727  Bool_t found(true);
728  if(!dynamic_cast<TH1 *>(testmgr.GetListOfHistograms()->FindObject("Test1D"))){
729  std::cout << "Not found: Test1D" << std::endl;
730  found = false;
731  }
732  if(!dynamic_cast<TH2 *>(testmgr.GetListOfHistograms()->FindObject("Test2D"))){
733  std::cout << "Not found: Test2D" << std::endl;
734  found = false;
735  }
736  if(!dynamic_cast<TH3 *>(testmgr.GetListOfHistograms()->FindObject("Test3D"))){
737  std::cout << "Not found: Test3D" << std::endl;
738  found = false;
739  }
740  if(!dynamic_cast<THnSparse *>(testmgr.GetListOfHistograms()->FindObject("TestNSparse"))){
741  std::cout << "Not found: TestNSparse" << std::endl;
742  found = false;
743  }
744  if(!dynamic_cast<TProfile *>(testmgr.GetListOfHistograms()->FindObject("TestProfile"))){
745  std::cout << "Not found: TestProfile" << std::endl;
746  found = false;
747  }
748 
749  return found ? 0 : 1;
750  }
751 
752  int THistManagerTestSuite::TestBuildGroupedHistograms(){
753  THistManager testmgr("testmgr");
754 
755  // Creating 3 groups
756  testmgr.CreateTH1("Group1/Test1", "Test Histogram 1 in group 1", 1, 0., 1.);
757  testmgr.CreateTH1("Group1/Test2", "Test Histogram 2 in group 1", 1, 0., 1.);
758  testmgr.CreateTH2("Group2/Test1", "Test Histogram 1 in group 2", 1, 0., 1., 1, 0., 1.);
759  testmgr.CreateTH2("Group2/Test2", "Test Histogram 2 in group 2", 1, 0., 1., 1, 0., 1.);
760  testmgr.CreateTH3("Group3/Test1", "Test Histogram 1 in group 3", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
761  testmgr.CreateTH3("Group3/Test2", "Test Histogram 2 in group 3", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
762  // Creating histogram in a subgroup
763  testmgr.CreateTProfile("Group4/Subgroup1/Test1", "Test histogram for subgroup handling", 1, 0., 1);
764 
765  // Evalutate test
766  // Tell user why test has failed
767  Bool_t found(true);
768  THashList *currentdir(nullptr);
769  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group1")))){
770  std::cout << "Not found: Group1" << std::endl;
771  found = false;
772  } else {
773  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test1"))){
774  std::cout << "Not found in Group1: Test1" << std::endl;
775  found = false;
776  }
777  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test2"))){
778  std::cout << "Not found in Group1: Test2" << std::endl;
779  found = false;
780  }
781  }
782  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group2")))){
783  std::cout << "Not found: Group2" << std::endl;
784  found = false;
785  } else {
786  if(!dynamic_cast<TH2 *>(currentdir->FindObject("Test1"))){
787  std::cout << "Not found in Group2: Test1" << std::endl;
788  found = false;
789  }
790  if(!dynamic_cast<TH2 *>(currentdir->FindObject("Test2"))){
791  std::cout << "Not found in Group2: Test2" << std::endl;
792  found = false;
793  }
794  }
795  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group3")))){
796  std::cout << "Not found: Group3" << std::endl;
797  found = false;
798  } else {
799  if(!static_cast<TH3 *>(currentdir->FindObject("Test1"))){
800  std::cout << "Not found in Group3: Test1" << std::endl;
801  found = false;
802  }
803  if(!static_cast<TH3 *>(currentdir->FindObject("Test2"))){
804  std::cout << "Not found in Group3: Test2" << std::endl;
805  found = false;
806  }
807  }
808  if(!(currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group4")))){
809  std::cout << "Not found: Group4" << std::endl;
810  found = false;
811  } else {
812  if(!(currentdir = dynamic_cast<THashList *>(currentdir->FindObject("Subgroup1")))){
813  std::cout << "Not found in Group4: Subgroup1" << std::endl;
814  found = false;
815  } else {
816  if(!dynamic_cast<TH1 *>(currentdir->FindObject("Test1"))){
817  std::cout << "Not found in Subgroup1: Test1" << std::endl;
818  found = false;
819  }
820  }
821  }
822  return found ? 0 : 1;
823  }
824 
825  int THistManagerTestSuite::TestFillSimpleHistograms(){
826  THistManager testmgr("testmgr");
827 
828  testmgr.CreateTH1("Test1", "Test fill 1D histogram", 1, 0., 1.);
829  testmgr.CreateTH2("Test2", "Test fill 2D histogram", 1, 0., 1., 1, 0., 1.);
830  testmgr.CreateTH3("Test3", "Test fill 3D histogram", 1, 0., 1., 1, 0., 1., 1, 0., 1.);
831  int nbins[4] = {1,1,1,1}; double min[4] = {0.,0.,0.,0.}, max[4] = {1.,1.,1.,1.};
832  testmgr.CreateTHnSparse("TestN", "Test Fill THnSparse", 4, nbins, min, max);
833  testmgr.CreateTProfile("TestProfile", "Test fill Profile histogram", 1, 0., 1.);
834 
835  double point[4] = {0.5, 0.5, 0.5, 0.5};
836  for(int i = 0; i < 100; i++){
837  testmgr.FillTH1("Test1", 0.5);
838  testmgr.FillTH2("Test2", 0.5, 0.5);
839  testmgr.FillTH3("Test3", 0.5, 0.5, 0.5);
840  testmgr.FillProfile("TestProfile", 0.5, 1.);
841  testmgr.FillTHnSparse("TestN", point);
842  }
843 
844  // Evalutate test
845  // tell user why test has failed
846  bool success(true);
847 
848  TH1 *test1 = dynamic_cast<TH1 *>(testmgr.GetListOfHistograms()->FindObject("Test1"));
849  if(test1){
850  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
851  std::cout << "Test1: Mismatch in values, expected 100, found " << test1->GetBinContent(1) << std::endl;
852  success = false;
853  }
854  } else {
855  std::cout << "Not found: Test1" << std::endl;
856  success = false;
857  }
858 
859  TH2 *test2 = dynamic_cast<TH2 *>(testmgr.GetListOfHistograms()->FindObject("Test2"));
860  if(test2){
861  if(TMath::Abs(test2->GetBinContent(1, 1) - 100) > DBL_EPSILON){
862  std::cout << "Test2: Mismatch in values, expected 100, found " << test2->GetBinContent(1,1) << std::endl;
863  success = false;
864  }
865  } else {
866  std::cout << "Not found: Test2" << std::endl;
867  success = false;
868  }
869 
870  TH3 *test3 = dynamic_cast<TH3 *>(testmgr.GetListOfHistograms()->FindObject("Test3"));
871  if(test3){
872  if(TMath::Abs(test3->GetBinContent(1, 1, 1) - 100) > DBL_EPSILON){
873  std::cout << "Test3: Mismatch in values, expected 100, found " << test3->GetBinContent(1,1,1) << std::endl;
874  success = false;
875  }
876  } else {
877  std::cout << "Not found: Test3" << std::endl;
878  success = false;
879  }
880 
881  THnSparse *testN = dynamic_cast<THnSparse *>(testmgr.GetListOfHistograms()->FindObject("TestN"));
882  if(testN){
883  int index[4] = {1,1,1,1};
884  if(TMath::Abs(testN->GetBinContent(index) - 100) > DBL_EPSILON){
885  std::cout << "TestN: Mismatch in values, expected 100, found " << testN->GetBinContent(index) << std::endl;
886  success = false;
887  }
888  } else {
889  std::cout << "Not found: TestN" << std::endl;
890  success = false;
891  }
892 
893  TProfile *testProfile = dynamic_cast<TProfile *>(testmgr.GetListOfHistograms()->FindObject("TestProfile"));
894  if(testProfile){
895  if(TMath::Abs(testProfile->GetBinContent(1) - 1) > DBL_EPSILON){
896  std::cout << "TestProfile: Mismatch in values, expected 1, found " << testProfile->GetBinContent(1) << std::endl;
897  success = false;
898  }
899  } else {
900  std::cout << "Not found: TestProfile" << std::endl;
901  success = false;
902  }
903 
904  return success ? 0 : 1;
905  }
906 
907 
908  int THistManagerTestSuite::TestFillGroupedHistograms(){
909  THistManager testmgr("testmgr");
910 
911  // Creating 3 groups, 1 with 1D and 1 with 2D histograms, and a third with a subgroup and a TProfile
912  testmgr.CreateTH1("Group1/Test1", "Test 1 Group 1D", 1, 0., 1.);
913  testmgr.CreateTH1("Group1/Test2", "Test 2 Group 1D", 1, 0., 1.);
914  testmgr.CreateTH2("Group2/Test1", "Test 1 Group 2D", 1, 0., 1., 1, 0., 1.);
915  testmgr.CreateTH2("Group2/Test2", "Test 2 Group 2D", 1, 0., 1., 1, 0., 1.);
916  testmgr.CreateTProfile("Group3/Subgroup1/Test1", "Test 1 with subgroup", 1, 0., 1.);
917 
918  for(int i = 0; i < 100; i++){
919  testmgr.FillTH1("Group1/Test1", 0.5);
920  testmgr.FillTH1("Group1/Test2", 0.5);
921  testmgr.FillTH2("Group2/Test1", 0.5, 0.5);
922  testmgr.FillTH2("Group2/Test2", 0.5, 0.5);
923  testmgr.FillProfile("Group3/Subgroup1/Test1", 0.5, 1);
924  }
925 
926  // Evaluate test
927  bool success(true);
928 
929  THashList *currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group1"));
930  if(currentdir){
931  TH1 *test1 = dynamic_cast<TH1 *>(currentdir->FindObject("Test1"));
932  if(test1){
933  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
934  std::cout << "Group1/Test1: Value mismatch: expected 100, found " << test1->GetBinContent(1) << std::endl;
935  success = false;
936  }
937  } else {
938  std::cout << "Not found in Group1: Test1" << std::endl;
939  success = false;
940  }
941  test1 = dynamic_cast<TH1 *>(currentdir->FindObject("Test2"));
942  if(test1){
943  if(TMath::Abs(test1->GetBinContent(1) - 100) > DBL_EPSILON){
944  std::cout << "Group1/Test2: Value mismatch: expected 100, found " << test1->GetBinContent(1) << std::endl;
945  success = false;
946  }
947  } else {
948  std::cout << "Not found in Group1: Test2" << std::endl;
949  success = false;
950  }
951  } else {
952  std::cout << "Not found: Group1" << std::endl;
953  success = false;
954  }
955 
956  currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group2"));
957  if(currentdir){
958  TH2 *test2 = dynamic_cast<TH2 *>(currentdir->FindObject("Test1"));
959  if(test2){
960  if(TMath::Abs(test2->GetBinContent(1,1) - 100) > DBL_EPSILON){
961  std::cout << "Group2/Test1: Value mismatch: expected 100, found " << test2->GetBinContent(1,1) << std::endl;
962  success = false;
963  }
964  } else {
965  std::cout << "Not found in Group2: Test1" << std::endl;
966  success = false;
967  }
968  test2 = dynamic_cast<TH2 *>(currentdir->FindObject("Test2"));
969  if(test2){
970  if(TMath::Abs(test2->GetBinContent(1,1) - 100) > DBL_EPSILON){
971  std::cout << "Group2/Test2: Value mismatch: expected 100, found " << test2->GetBinContent(1,1) << std::endl;
972  success = false;
973  }
974  } else {
975  std::cout << "Not found in Group2: Test2" << std::endl;
976  success = false;
977  }
978  } else {
979  std::cout << "Not found: Group2" << std::endl;
980  success = false;
981  }
982 
983  currentdir = dynamic_cast<THashList *>(testmgr.GetListOfHistograms()->FindObject("Group3"));
984  if(currentdir){
985  currentdir = dynamic_cast<THashList *>(currentdir->FindObject("Subgroup1"));
986  if(currentdir){
987  TProfile *testprofile = dynamic_cast<TProfile *>(currentdir->FindObject("Test1"));
988  if(testprofile){
989  if(TMath::Abs(testprofile->GetBinContent(1) - 1) > DBL_EPSILON){
990  std::cout << "Group3/Subgroup1/Test1: Value mismatch: expected 1, found " << testprofile->GetBinContent(1) << std::endl;
991  success = false;
992  }
993  } else {
994  std::cout << "Not found in Group3/Subgroup1: Test1" << std::endl;
995  success = false;
996  }
997  } else {
998  std::cout << "Not found in Group3: Subgroup1" << std::endl;
999  success = false;
1000  }
1001  } else {
1002  std::cout << "Not found: Group3" << std::endl;
1003  success = false;
1004  }
1005  return success ? 0 : 1;
1006  }
1007 
1008  int TestRunAll(){
1009  int testresult(0);
1010  THistManagerTestSuite testsuite;
1011 
1012  std::cout << "Running test: Build simple" << std::endl;
1013  testresult += testsuite.TestBuildSimpleHistograms();
1014  std::cout << "Result after test: " << testresult << std::endl;
1015 
1016  std::cout << "Running test: Build grouped" << std::endl;
1017  testresult += testsuite.TestBuildGroupedHistograms();
1018  std::cout << "Result after test: " << testresult << std::endl;
1019 
1020  std::cout << "Running test: Fill Simple" << std::endl;
1021  testresult += testsuite.TestFillSimpleHistograms();
1022  std::cout << "Result after test: " << testresult << std::endl;
1023 
1024  std::cout << "Running test: Fill Grouped" << std::endl;
1025  testresult += testsuite.TestFillGroupedHistograms();
1026  std::cout << "Result after test: " << testresult << std::endl;
1027 
1028  return testresult;
1029  }
1030 
1032  THistManagerTestSuite testsuite;
1033  return testsuite.TestBuildSimpleHistograms();
1034  }
1035 
1037  THistManagerTestSuite testsuite;
1038  return testsuite.TestBuildGroupedHistograms();
1039  }
1040 
1042  THistManagerTestSuite testsuite;
1043  return testsuite.TestFillSimpleHistograms();
1044  }
1045 
1047  THistManagerTestSuite testsuite;
1048  return testsuite.TestFillGroupedHistograms();
1049  }
1050 }
Int_t fNext
Next position in the histmanager.
Definition: THistManager.h:226
iterator & operator=(const iterator &rhs)
THashList * CreateHistoGroup(const char *groupname)
Create a new group of histograms within a parent group.
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:710
Int_t nbinsy
THistManager()
Default constructor.
THashList * FindGroup(const char *dirname) const
Find histogram group.
void FillTH2(const char *hname, double x, double y, double weight=1., Option_t *opt="")
Fill a 2D histogram within the container.
void FillTH3(const char *hname, double x, double y, double z, double weight=1., Option_t *opt="")
Fill a 3D histogram within the container.
virtual void CreateBinEdges(TArrayD &binedges) const =0
const THistManager * fkArray
Underlying histmanager to iterate over.
Definition: THistManager.h:224
void SetObject(TObject *const o, const char *group="/")
Set a new group into the container into the parent group.
Interface for binnings used by the histogram handler.
Definition: TBinning.h:21
THMIDirection_t fDirection
Direction of the iterator.
Definition: THistManager.h:227
Bool_t operator!=(const iterator &aIter) const
Comparison operator for unequalness.
void FillTHnSparse(const char *name, const double *x, double weight=1., Option_t *opt="")
~THistManager()
Destructor.
TH2 * CreateTH2(const char *name, const char *title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, Option_t *opt="")
Create a new TH2 within the container.
TObject * FindObject(const char *name) const
Find an object inside the container.
int Int_t
Definition: External.C:63
void CreateTProfile(const char *name, const char *title, int nbinsX, double xmin, double xmax, Option_t *opt="")
Create a new TProfile within the container.
THashList * GetListOfHistograms() const
Get the list of histograms.
Definition: THistManager.h:658
Definition: External.C:252
TH1 * CreateTH1(const char *name, const char *title, int nbins, double xmin, double xmax, Option_t *opt="")
Create a new TH1 within the container.
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:128
stl-iterator for the histogram manager
Definition: THistManager.h:115
void FillTH1(const char *hname, double x, double weight=1., Option_t *opt="")
Fill a 1D histogram within the container.
iterator & operator++()
Prefix increment operator.
Int_t nbinsx
TString basename(const TString &path) const
Extracting the basename from a given histogram path.
Collection of tests for the THistManager.
Definition: THistManager.h:753
iterator & operator--()
Prefix decrement operator.
TString histname(const TString &path) const
Extracting the histogram name from a given histogram path.
Definition: External.C:220
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
TObject * operator*() const
Dereferncing operator.
const Double_t ymin
Definition: AddTaskCFDStar.C:6
Container class for histograms.
Definition: THistManager.h:99
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="")
Create a new THnSparse within the container.
THashList * fHistos
List of histograms.
Definition: THistManager.h:709
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="")
Create a new TH2 within the container.
Int_t fCurrentPos
Current position of the iterator in the histmanager.
Definition: THistManager.h:225