AliPhysics  e34b7ac (e34b7ac)
 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 <cstring>
16 #include <sstream>
17 #include <string>
18 #include <exception>
19 #include <vector>
20 #include <TArrayD.h>
21 #include <TAxis.h>
22 #include <TH1.h>
23 #include <TH2.h>
24 #include <TH3.h>
25 #include <THnSparse.h>
26 #include <THashList.h>
27 #include <TObjArray.h>
28 #include <TObjString.h>
29 #include <TProfile.h>
30 #include <TString.h>
31 
32 #include "TBinning.h"
33 #include "THistManager.h"
34 
38 
40  TNamed(),
41  fHistos(NULL),
42  fIsOwner(true)
43 {
44 }
45 
46 THistManager::THistManager(const char *name):
47  TNamed(name, Form("Histogram container %s", name)),
48  fHistos(NULL),
49  fIsOwner(true)
50 {
51  fHistos = new THashList();
52  fHistos->SetName(Form("histos%s", name));
53  fHistos->SetOwner();
54 }
55 
57  if(fHistos && fIsOwner) delete fHistos;
58 }
59 
60 THashList* THistManager::CreateHistoGroup(const char *groupname, const char *parent) {
61  THashList *parentgroup = FindGroup(parent);
62  if(!parentgroup){
63  Fatal("THistManager::CreateHistoGroup", "Parent group %s does not exist", parent);
64  return 0;
65  }
66  THashList *childgroup = new THashList();
67  childgroup->SetName(groupname);
68  parentgroup->Add(childgroup);
69  return childgroup;
70 }
71 
72 TH1* THistManager::CreateTH1(const char *name, const char *title, int nbins, double xmin, double xmax, Option_t *opt){
73  TString dirname(basename(name)), hname(histname(name));
74  THashList *parent(FindGroup(dirname.Data()));
75  if(!parent){
76  Fatal("THistManager::CreateTH1", "Parent %s does not exist", dirname.Data());
77  return 0 ;
78  }
79  if(parent->FindObject(hname.Data())){
80  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
81  return 0;
82  }
83  TH1* h = new TH1D(hname.Data(), title, nbins, xmin, xmax);
84  TString optionstring(opt);
85  optionstring.ToLower();
86  if(optionstring.Contains("s"))
87  h->Sumw2();
88  parent->Add(h);
89  return h;
90 }
91 
92 TH1* THistManager::CreateTH1(const char *name, const char *title, int nbins, const double *xbins, Option_t *opt){
93  TString dirname(basename(name)), hname(histname(name));
94  THashList *parent(FindGroup(dirname));
95  if(!parent){
96  Fatal("THistManager::CreateTH1", "Parent %s does not exist", dirname.Data());
97  return 0;
98  }
99  if(parent->FindObject(hname.Data())){
100  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
101  return 0;
102  }
103  TH1* h = new TH1D(hname.Data(), title, nbins, xbins);
104  TString optionstring(opt);
105  optionstring.ToLower();
106  if(optionstring.Contains("s"))
107  h->Sumw2();
108  parent->Add(h);
109  return h;
110 }
111 
112 TH1* THistManager::CreateTH1(const char *name, const char *title, const TArrayD &xbins, Option_t *opt){
113  TString dirname(basename(name)), hname(histname(name));
114  THashList *parent(FindGroup(dirname));
115  if(!parent){
116  Fatal("THistManager::CreateTH1", "Parent %s does not exist", dirname.Data());
117  return 0;
118  }
119  if(parent->FindObject(hname.Data())){
120  Fatal("THistManager::CreateTH1", "Object %s already exists in group %s", hname.Data(), dirname.Data());
121  return 0;
122  }
123  TH1* h = new TH1D(hname.Data(), title, xbins.GetSize()-1, xbins.GetArray());
124  TString optionstring(opt);
125  optionstring.ToLower();
126  if(optionstring.Contains("s"))
127  h->Sumw2();
128  parent->Add(h);
129  return h;
130 }
131 
132 TH1* THistManager::CreateTH1(const char *name, const char *title, const TBinning &xbin, Option_t *opt){
133  TArrayD myxbins;
134  try{
135  xbin.CreateBinEdges(myxbins);
136  } catch(std::exception &e){
137  Fatal("THistManager::CreateTH1", "Exception raised: %s", e.what());
138  }
139  return CreateTH1(name, title, myxbins, opt);
140 }
141 
142 TH2* THistManager::CreateTH2(const char *name, const char *title, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, Option_t *opt){
143  TString dirname(basename(name)), hname(histname(name));
144  THashList *parent(FindGroup(dirname.Data()));
145  if(!parent){
146  Fatal("THistManager::CreateTH2", "Parent %s does not exist", dirname.Data());
147  return 0;
148  }
149  if(parent->FindObject(hname.Data())){
150  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
151  return 0;
152  }
153  TH2* h = new TH2D(hname.Data(), title, nbinsx, xmin, xmax, nbinsy, ymin, ymax);
154  TString optionstring(opt);
155  optionstring.ToLower();
156  if(optionstring.Contains("s"))
157  h->Sumw2();
158  parent->Add(h);
159  return h;
160 }
161 
162 TH2* THistManager::CreateTH2(const char *name, const char *title, int nbinsx, const double *xbins, int nbinsy, const double *ybins, Option_t *opt){
163  TString dirname(basename(name)), hname(histname(name));
164  THashList *parent(FindGroup(dirname.Data()));
165  if(!parent){
166  Fatal("THistManager::CreateTH2", "Parent %s does not exist", dirname.Data());
167  return 0;
168  }
169  if(parent->FindObject(hname.Data())){
170  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
171  return 0;
172  }
173  TH2* h = new TH2D(hname.Data(), title, nbinsx, xbins, nbinsy, ybins);
174  TString optionstring(opt);
175  optionstring.ToLower();
176  if(optionstring.Contains("s"))
177  h->Sumw2();
178  parent->Add(h);
179  return h;
180 }
181 
182 TH2* THistManager::CreateTH2(const char *name, const char *title, const TArrayD &xbins, const TArrayD &ybins, Option_t *opt){
183  TString dirname(basename(name)), hname(histname(name));
184  THashList *parent(FindGroup(dirname.Data()));
185  if(!parent){
186  Fatal("THistManager::CreateTH2", "Parent %s does not exist", dirname.Data());
187  return 0;
188  }
189  if(parent->FindObject(hname.Data())){
190  Fatal("THistManager::CreateTH2", "Object %s already exists in group %s", hname.Data(), dirname.Data());
191  return 0;
192  }
193  TH2* h = new TH2D(hname.Data(), title, xbins.GetSize() - 1, xbins.GetArray(), ybins.GetSize() - 1, ybins.GetArray());
194  TString optionstring(opt);
195  optionstring.ToLower();
196  if(optionstring.Contains("s"))
197  h->Sumw2();
198  parent->Add(h);
199  return h;
200 }
201 
202 TH2* THistManager::CreateTH2(const char *name, const char *title, const TBinning &xbins, const TBinning &ybins, Option_t *opt){
203  TArrayD myxbins, myybins;
204  try{
205  xbins.CreateBinEdges(myxbins);
206  } catch (std::exception &e) {
207  Fatal("THistManager::CreateTH2 (x-dir)", "Exception raised: %s", e.what());
208  }
209  try{
210  ybins.CreateBinEdges(myybins);
211  } catch (std::exception &e) {
212  Fatal("THistManager::CreateTH2 (y-dir)", "Exception raised: %s", e.what());
213  }
214 
215  return CreateTH2(name, title, myxbins, myybins, opt);
216 }
217 
218 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) {
219  TString dirname(basename(name)), hname(histname(name));
220  THashList *parent(FindGroup(dirname.Data()));
221  if(!parent){
222  Fatal("THistManager::CreateTH3", "Parent %s does not exist", dirname.Data());
223  return 0;
224  }
225  if(parent->FindObject(hname.Data())){
226  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
227  return 0;
228  }
229  TH3* h = new TH3D(hname.Data(), title, nbinsx, xmin, xmax, nbinsy, ymin, ymax, nbinsz, zmin, zmax);
230  TString optionstring(opt);
231  optionstring.ToLower();
232  if(optionstring.Contains("s"))
233  h->Sumw2();
234  parent->Add(h);
235  return h;
236 }
237 
238 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) {
239  TString dirname(basename(name)), hname(histname(name));
240  THashList *parent(FindGroup(dirname.Data()));
241  if(!parent){
242  Fatal("THistManager::CreateTH3", "Parent %s does not exist", dirname.Data());
243  return 0;
244  }
245  if(parent->FindObject(hname.Data())){
246  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
247  return 0;
248  }
249  TH3* h = new TH3D(hname.Data(), title, nbinsx, xbins, nbinsy, ybins, nbinsz, zbins);
250  TString optionstring(opt);
251  optionstring.ToLower();
252  if(optionstring.Contains("s"))
253  h->Sumw2();
254  parent->Add(h);
255  return h;
256 }
257 
258 TH3* THistManager::CreateTH3(const char* name, const char* title, const TArrayD& xbins, const TArrayD& ybins, const TArrayD& zbins, Option_t *opt) {
259  TString dirname(basename(name)), hname(histname(name));
260  THashList *parent(FindGroup(dirname.Data()));
261  if(!parent){
262  Fatal("THistManager::CreateTH3", "Parent %s does not exist", dirname.Data());
263  return 0;
264  }
265  if(parent->FindObject(hname.Data())){
266  Fatal("THistManager::CreateTH3", "Object %s already exists in group %s", hname.Data(), dirname.Data());
267  return 0;
268  }
269  TH3* h = new TH3D(hname.Data(), title, xbins.GetSize()-1, xbins.GetArray(), ybins.GetSize()-1, ybins.GetArray(), zbins.GetSize()-1, zbins.GetArray());
270  TString optionstring(opt);
271  optionstring.ToLower();
272  if(optionstring.Contains("s"))
273  h->Sumw2();
274  parent->Add(h);
275  return h;
276 }
277 
278 TH3* THistManager::CreateTH3(const char *name, const char *title, const TBinning &xbins, const TBinning &ybins, const TBinning &zbins, Option_t *opt){
279  TArrayD myxbins, myybins, myzbins;
280  try{
281  xbins.CreateBinEdges(myxbins);
282  } catch (std::exception &e) {
283  Fatal("THistManager::CreateTH2 (x-dir)", "Exception raised: %s", e.what());
284  }
285  try{
286  ybins.CreateBinEdges(myybins);
287  } catch (std::exception &e) {
288  Fatal("THistManager::CreateTH2 (y-dir)", "Exception raised: %s", e.what());
289  }
290  try{
291  zbins.CreateBinEdges(myzbins);
292  } catch (std::exception &e) {
293  Fatal("THistManager::CreateTH2 (z-dir)", "Exception raised: %s", e.what());
294  }
295 
296  return CreateTH3(name, title, myxbins, myybins, myzbins);
297 }
298 
299 THnSparse* THistManager::CreateTHnSparse(const char *name, const char *title, int ndim, const int *nbins, const double *min, const double *max, Option_t *opt) {
300  TString dirname(basename(name)), hname(histname(name));
301  THashList *parent(FindGroup(dirname.Data()));
302  if(!parent){
303  Fatal("THistManager::CreateTHnSparse", "Parent %s does not exist", dirname.Data());
304  return 0;
305  }
306  if(parent->FindObject(hname.Data())){
307  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
308  return 0;
309  }
310  THnSparse* h = new THnSparseD(hname.Data(), title, ndim, nbins, min, max);
311  TString optionstring(opt);
312  optionstring.ToLower();
313  if(optionstring.Contains("s"))
314  h->Sumw2();
315  parent->Add(h);
316  return h;
317 }
318 
319 THnSparse* THistManager::CreateTHnSparse(const char *name, const char *title, int ndim, const TAxis **axes, Option_t *opt) {
320  TString dirname(basename(name)), hname(histname(name));
321  THashList *parent(FindGroup(dirname.Data()));
322  if(!parent){
323  Fatal("THistManager::CreateTHnSparse", "Parent %s does not exist", dirname.Data());
324  return 0;
325  }
326  if(parent->FindObject(hname)){
327  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
328  return 0;
329  }
330  TArrayD xmin(ndim), xmax(ndim);
331  TArrayI nbins(ndim);
332  for(int idim = 0; idim < ndim; ++idim){
333  const TAxis &myaxis = *(axes[idim]);
334  nbins[idim] = myaxis.GetNbins();
335  xmin[idim] = myaxis.GetXmin();
336  xmax[idim] = myaxis.GetXmax();
337  }
338  THnSparseD *hsparse = new THnSparseD(hname.Data(), title, ndim, nbins.GetArray(), xmin.GetArray(), xmax.GetArray());
339  for(int id = 0; id < ndim; ++id)
340  *(hsparse->GetAxis(id)) = *(axes[id]);
341  TString optionstring(opt);
342  optionstring.ToLower();
343  if(optionstring.Contains("s"))
344  hsparse->Sumw2();
345  parent->Add(hsparse);
346  return hsparse;
347 }
348 
349 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, double xmin, double xmax, Option_t *opt) {
350  TString dirname(basename(name)), hname(histname(name));
351  THashList *parent(FindGroup(dirname));
352  if(!parent)
353  Fatal("THistManager::CreateTProfile", "Parent %s does not exist", dirname.Data());
354  if(parent->FindObject(hname.Data()))
355  Fatal("THistManager::CreateTProfile", "Object %s already exists in group %s", hname.Data(), dirname.Data());
356  TProfile *hist = new TProfile(hname.Data(), title, nbinsX, xmin, xmax, opt);
357  parent->Add(hist);
358 }
359 
360 void THistManager::CreateTProfile(const char* name, const char* title, int nbinsX, const double* xbins, Option_t *opt) {
361  TString dirname(basename(name)), hname(histname(name));
362  THashList *parent(FindGroup(dirname));
363  if(!parent)
364  Fatal("THistManager::CreateTHnSparse", "Parent %s does not exist", dirname.Data());
365  if(parent->FindObject(hname.Data()))
366  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
367  TProfile *hist = new TProfile(hname.Data(), title, nbinsX, xbins, opt);
368  parent->Add(hist);
369 }
370 
371 void THistManager::CreateTProfile(const char* name, const char* title, const TArrayD& xbins, Option_t *opt){
372  TString dirname(basename(name)), hname(histname(name));
373  THashList *parent(FindGroup(dirname));
374  if(!parent)
375  Fatal("THistManager::CreateTHnSparse", "Parent %s does not exist", dirname.Data());
376  if(parent->FindObject(hname.Data()))
377  Fatal("THistManager::CreateTHnSparse", "Object %s already exists in group %s", hname.Data(), dirname.Data());
378  TProfile *hist = new TProfile(hname.Data(), title, xbins.GetSize()-1, xbins.GetArray(), opt);
379  parent->Add(hist);
380 }
381 
382 void THistManager::CreateTProfile(const char *name, const char *title, const TBinning &xbins, Option_t *opt){
383  TArrayD myxbins;
384  try{
385  xbins.CreateBinEdges(myxbins);
386  } catch (std::exception &e){
387  Fatal("THistManager::CreateProfile", "Exception raised: %s", e.what());
388  }
389  CreateTProfile(name, title, myxbins, opt);
390 }
391 
392 void THistManager::SetObject(TObject * const o, const char *group) {
393  THashList *parent(FindGroup(group));
394  if(!parent){
395  Fatal("THistManager::SetObject", "Parent %s does not exist", strcmp(group, "/") ? group : "");
396  return;
397  }
398  if(parent->FindObject(o->GetName())){
399  Fatal("THistManager::SetObject", "Parent %s does not exist", strcmp(group, "/") ? group : "");
400  return;
401  }
402  if(!(dynamic_cast<THnBase *>(o) || dynamic_cast<TH1 *>(o))){
403  Fatal("THistManager::SetObject", "Object %s is not of a histogram type",o->GetName());
404  return;
405  }
406  fHistos->Add(o);
407 }
408 
409 void THistManager::FillTH1(const char *name, double x, double weight, Option_t *opt) {
410  TString dirname(basename(name)), hname(histname(name));
411  THashList *parent(FindGroup(dirname.Data()));
412  if(!parent){
413  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
414  return;
415  }
416  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname.Data()));
417  if(!hist){
418  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
419  return;
420  }
421  TString optionstring(opt);
422  if(optionstring.Contains("w")){
423  // use bin width as weight
424  Int_t bin = hist->GetXaxis()->FindBin(x);
425  // check if not overflow or underflow bin
426  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
427  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
428  }
429  hist->Fill(x, weight);
430 }
431 
432 void THistManager::FillTH1(const char *name, const char *label, double weight, Option_t *opt) {
433  TString dirname(basename(name)), hname(histname(name));
434  THashList *parent(FindGroup(dirname.Data()));
435  if(!parent){
436  Fatal("THistManager::FillTH1", "Parnt group %s does not exist", dirname.Data());
437  return;
438  }
439  TH1 *hist = dynamic_cast<TH1 *>(parent->FindObject(hname.Data()));
440  if(!hist){
441  Fatal("THistManager::FillTH1", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
442  return;
443  }
444  TString optionstring(opt);
445  if(optionstring.Contains("w")){
446  // use bin width as weight
447  // get bin for label
448  Int_t bin = hist->GetXaxis()->FindBin(label);
449  // check if not overflow or underflow bin
450  if(bin != 0 && bin != hist->GetXaxis()->GetNbins())
451  weight = 1./hist->GetXaxis()->GetBinWidth(bin);
452  }
453  hist->Fill(label, weight);
454 }
455 
456 void THistManager::FillTH2(const char *name, double x, double y, double weight, Option_t *opt) {
457  TString dirname(basename(name)), hname(histname(name));
458  THashList *parent(FindGroup(dirname.Data()));
459  if(!parent){
460  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
461  return;
462  }
463  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname.Data()));
464  if(!hist){
465  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
466  return;
467  }
468  TString optstring(opt);
469  Double_t myweight = optstring.Contains("w") ? 1. : weight;
470  if(optstring.Contains("wx")){
471  Int_t binx = hist->GetXaxis()->FindBin(x);
472  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
473  }
474  if(optstring.Contains("wy")){
475  Int_t biny = hist->GetYaxis()->FindBin(y);
476  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
477  }
478  hist->Fill(x, y, myweight);
479 }
480 
481 void THistManager::FillTH2(const char *name, double *point, double weight, Option_t *opt) {
482  TString dirname(basename(name)), hname(histname(name));
483  THashList *parent(FindGroup(dirname.Data()));
484  if(!parent){
485  Fatal("THistManager::FillTH2", "Parent group %s does not exist", dirname.Data());
486  return;
487  }
488  TH2 *hist = dynamic_cast<TH2 *>(parent->FindObject(hname.Data()));
489  if(!hist){
490  Fatal("THistManager::FillTH2", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
491  return;
492  }
493  TString optstring(opt);
494  Double_t myweight = optstring.Contains("w") ? 1. : weight;
495  if(optstring.Contains("wx")){
496  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
497  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
498  }
499  if(optstring.Contains("wy")){
500  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
501  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
502  }
503  hist->Fill(point[0], point[1], weight);
504 }
505 
506 void THistManager::FillTH3(const char* name, double x, double y, double z, double weight, Option_t *opt) {
507  TString dirname(basename(name)), hname(histname(name));
508  THashList *parent(FindGroup(dirname.Data()));
509  if(!parent){
510  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
511  return;
512  }
513  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname.Data()));
514  if(!hist){
515  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
516  return;
517  }
518  TString optstring(opt);
519  Double_t myweight = optstring.Contains("w") ? 1. : weight;
520  if(optstring.Contains("wx")){
521  Int_t binx = hist->GetXaxis()->FindBin(x);
522  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
523  }
524  if(optstring.Contains("wy")){
525  Int_t biny = hist->GetYaxis()->FindBin(y);
526  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
527  }
528  if(optstring.Contains("wz")){
529  Int_t binz = hist->GetZaxis()->FindBin(z);
530  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
531  }
532  hist->Fill(x, y, z, weight);
533 }
534 
535 void THistManager::FillTH3(const char* name, const double* point, double weight, Option_t *opt) {
536  TString dirname(basename(name)), hname(histname(name));
537  THashList *parent(FindGroup(dirname.Data()));
538  if(!parent){
539  Fatal("THistManager::FillTH3", "Parent group %s does not exist", dirname.Data());
540  return;
541  }
542  TH3 *hist = dynamic_cast<TH3 *>(parent->FindObject(hname.Data()));
543  if(!hist){
544  Fatal("THistManager::FillTH3", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
545  return;
546  }
547  TString optstring(opt);
548  Double_t myweight = optstring.Contains("w") ? 1. : weight;
549  if(optstring.Contains("wx")){
550  Int_t binx = hist->GetXaxis()->FindBin(point[0]);
551  if(binx != 0 && binx != hist->GetXaxis()->GetNbins()) myweight *= 1./hist->GetXaxis()->GetBinWidth(binx);
552  }
553  if(optstring.Contains("wy")){
554  Int_t biny = hist->GetYaxis()->FindBin(point[1]);
555  if(biny != 0 && biny != hist->GetYaxis()->GetNbins()) myweight *= 1./hist->GetYaxis()->GetBinWidth(biny);
556  }
557  if(optstring.Contains("wz")){
558  Int_t binz = hist->GetZaxis()->FindBin(point[2]);
559  if(binz != 0 && binz != hist->GetZaxis()->GetNbins()) myweight *= 1./hist->GetZaxis()->GetBinWidth(binz);
560  }
561  hist->Fill(point[0], point[1], point[2], weight);
562 }
563 
564 void THistManager::FillTHnSparse(const char *name, const double *x, double weight, Option_t *opt) {
565  TString dirname(basename(name)), hname(histname(name));
566  THashList *parent(FindGroup(dirname.Data()));
567  if(!parent){
568  Fatal("THistManager::FillTHnSparse", "Parent group %s does not exist", dirname.Data());
569  return;
570  }
571  THnSparseD *hist = dynamic_cast<THnSparseD *>(parent->FindObject(hname.Data()));
572  if(!hist){
573  Fatal("THistManager::FillTHnSparse", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
574  return;
575  }
576  TString optstring(opt);
577  Double_t myweight = optstring.Contains("w") ? 1. : weight;
578  for(Int_t iaxis = 0; iaxis < hist->GetNdimensions(); iaxis++){
579  std::stringstream weighthandler;
580  weighthandler << "w" << iaxis;
581  if(optstring.Contains(weighthandler.str().c_str())){
582  Int_t bin = hist->GetAxis(iaxis)->FindBin(x[iaxis]);
583  if(bin != 0 && bin != hist->GetAxis(iaxis)->GetNbins()) myweight *= hist->GetAxis(iaxis)->GetBinWidth(bin);
584  }
585  }
586 
587  hist->Fill(x, weight);
588 }
589 
590 void THistManager::FillProfile(const char* name, double x, double y, double weight){
591  TString dirname(basename(name)), hname(histname(name));
592  THashList *parent(FindGroup(dirname.Data()));
593  if(!parent)
594  Fatal("THistManager::FillTProfile", "Parent group %s does not exist", dirname.Data());
595  TProfile *hist = dynamic_cast<TProfile *>(parent->FindObject(hname.Data()));
596  if(!hist)
597  Fatal("THistManager::FillTProfile", "Histogram %s not found in parent group %s", hname.Data(), dirname.Data());
598  hist->Fill(x, y, weight);
599 }
600 
601 TObject *THistManager::FindObject(const char *name) const {
602  TString dirname(basename(name)), hname(histname(name));
603  THashList *parent(FindGroup(dirname.Data()));
604  if(!parent) return NULL;
605  return parent->FindObject(hname);
606 }
607 
609  TString dirname(basename(obj->GetName())), hname(histname(obj->GetName()));
610  THashList *parent(FindGroup(dirname.Data()));
611  if(!parent) return NULL;
612  return parent->FindObject(hname);
613 }
614 
615 THashList *THistManager::FindGroup(const char *dirname) const {
616  if(!strlen(dirname) || !strcmp(dirname, "/")) return fHistos;
617  std::vector<std::string> tokens;
618  TokenizeFilename(dirname, "/", tokens);
619  THashList *currentdir(fHistos);
620  for(std::vector<std::string>::iterator it = tokens.begin(); it != tokens.end(); ++it){
621  currentdir = dynamic_cast<THashList *>(currentdir->FindObject(it->c_str()));
622  if(!currentdir) break;
623  }
624  return currentdir;
625 }
626 
627 void THistManager::TokenizeFilename(const char *name, const char *delim, std::vector<std::string> &listoftokens) const {
628  TString s(name);
629  TObjArray *arr = s.Tokenize(delim);
630  TObjString *ostr(NULL);
631  TIter toks(arr);
632  while((ostr = dynamic_cast<TObjString *>(toks()))){
633  listoftokens.push_back(std::string(ostr->String().Data()));
634  }
635  delete arr;
636 }
637 
638 TString THistManager::basename(const char *path) const {
639  TString s(path);
640  int index = s.Last('/');
641  if(index < 0) return ""; // no directory structure
642  return TString(s(0, index)).Data();
643 }
644 
645 TString THistManager::histname(const char *path) const {
646  TString s(path);
647  int index = s.Last('/');
648  if(index < 0) return path; // no directory structure
649  return TString(s(index+1, s.Length() - (index+1))).Data();
650 }
651 
657 
659  fkArray(hmgr),
660  fCurrentPos(),
661  fNext(),
662  fDirection(dir)
663 {}
664 
666  fkArray(ref.fkArray),
667  fCurrentPos(ref.fCurrentPos),
668  fNext(ref.fNext),
669  fDirection(ref.fDirection)
670 {}
671 
673  if(this != &ref){
674  fkArray = ref.fkArray;
675  fCurrentPos = ref.fCurrentPos;
676  fNext = ref.fNext;
677  fDirection = ref.fDirection;
678  }
679  return *this;
680 }
681 
683  return fCurrentPos == other.fCurrentPos;
684 }
685 
687  if(fDirection == kTHMIforward)
688  fCurrentPos++;
689  else
690  fCurrentPos--;
691  return *this;
692 }
693 
695  iterator tmp(*this);
696  operator++();
697  return tmp;
698 }
699 
701  if(fDirection == kTHMIforward)
702  fCurrentPos--;
703  else
704  fCurrentPos++;
705  return *this;
706 };
707 
709  iterator tmp(*this);
710  operator--();
711  return tmp;
712 };
713 
715  if(fCurrentPos >=0 && fCurrentPos < fkArray->GetListOfHistograms()->GetEntries())
716  return fkArray->GetListOfHistograms()->At(fCurrentPos);
717  return NULL;
718 }
iterator & operator=(const iterator &rhs)
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:557
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
void TokenizeFilename(const char *name, const char *delim, std::vector< std::string > &listoftokens) const
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
TString histname(const char *path) const
THashList * CreateHistoGroup(const char *groupname, const char *parent="/")
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
TString basename(const char *path) const
stl-iterator for the histogram manager
Definition: THistManager.h:59
void FillTH1(const char *hname, double x, double weight=1., Option_t *opt="")
Bool_t Data(TH1F *h, Double_t *rangefit, Bool_t writefit, Double_t &sgn, Double_t &errsgn, Double_t &bkg, Double_t &errbkg, Double_t &sgnf, Double_t &errsgnf, Double_t &sigmafit, Int_t &status)
Int_t nbinsx
Definition: External.C:220
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
TObject * operator*() const
const Double_t ymin
Definition: AddTaskCFDStar.C:6
Container class for histograms for the high- charged particle analysis.
Definition: THistManager.h:43
const char Option_t
Definition: External.C:48
const Double_t zmax
const Int_t nbins
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:556
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="")