1 /*************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 // Generic Histogram container with support for groups and filling of groups by passing
21 // Jens Wiechula <Jens.Wiechula@cern.ch>
22 // Julian Book <Julian.Book@cern.ch>
31 #include <THnSparse.h>
33 #include <TProfile2D.h>
34 #include <TProfile3D.h>
35 #include <TCollection.h>
36 #include <THashList.h>
38 #include <TObjString.h>
39 #include <TObjArray.h>
48 #include <TVirtualPS.h>
51 #include "AliDielectronHelper.h"
52 #include "AliDielectronVarManager.h"
53 #include "AliDielectronHistos.h"
55 ClassImp(AliDielectronHistos)
58 AliDielectronHistos::AliDielectronHistos() :
60 TNamed("AliDielectronHistos","Dielectron Histogram Container"),
63 fUsedVars(new TBits(AliDielectronVarManager::kNMaxValues)),
64 fReservedWords(new TString)
67 // Default constructor
69 fHistoList.SetOwner(kTRUE);
70 fHistoList.SetName("Dielectron_Histos");
73 //_____________________________________________________________________________
74 AliDielectronHistos::AliDielectronHistos(const char* name, const char* title) :
79 fUsedVars(new TBits(AliDielectronVarManager::kNMaxValues)),
80 fReservedWords(new TString)
85 fHistoList.SetOwner(kTRUE);
86 fHistoList.SetName(name);
89 //_____________________________________________________________________________
90 AliDielectronHistos::~AliDielectronHistos()
96 if (fUsedVars) delete fUsedVars;
97 if (fList) fList->Clear();
98 delete fReservedWords;
101 //_____________________________________________________________________________
102 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
104 Int_t nbinsX, Double_t xmin, Double_t xmax,
105 UInt_t valTypeX, Bool_t logBinX, TString option,
109 // Default histogram creation 1D case
112 TVectorD *binLimX=0x0;
115 binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
117 binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
119 UserProfile(histClass,name,title,valTypeP,binLimX,valTypeX,option,valTypeW);
122 //_____________________________________________________________________________
123 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
125 Int_t nbinsX, Double_t xmin, Double_t xmax,
126 Int_t nbinsY, Double_t ymin, Double_t ymax,
127 UInt_t valTypeX, UInt_t valTypeY,
128 Bool_t logBinX, Bool_t logBinY, TString option,
132 // Default histogram creation 2D case
134 if (!IsHistogramOk(histClass,name)) return;
136 TVectorD *binLimX=0x0;
137 TVectorD *binLimY=0x0;
140 binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
142 binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
145 binLimY=AliDielectronHelper::MakeLogBinning(nbinsY, ymin, ymax);
147 binLimY=AliDielectronHelper::MakeLinBinning(nbinsY, ymin, ymax);
149 UserProfile(histClass,name,title,valTypeP,binLimX,binLimY,valTypeX,valTypeY,option,valTypeW);
153 //_____________________________________________________________________________
154 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
156 Int_t nbinsX, Double_t xmin, Double_t xmax,
157 Int_t nbinsY, Double_t ymin, Double_t ymax,
158 Int_t nbinsZ, Double_t zmin, Double_t zmax,
159 UInt_t valTypeX, UInt_t valTypeY, UInt_t valTypeZ,
160 Bool_t logBinX, Bool_t logBinY, Bool_t logBinZ, TString option,
164 // Default histogram creation 3D case
166 if (!IsHistogramOk(histClass,name)) return;
168 TVectorD *binLimX=0x0;
169 TVectorD *binLimY=0x0;
170 TVectorD *binLimZ=0x0;
173 binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
175 binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
179 binLimY=AliDielectronHelper::MakeLogBinning(nbinsY, ymin, ymax);
181 binLimY=AliDielectronHelper::MakeLinBinning(nbinsY, ymin, ymax);
185 binLimZ=AliDielectronHelper::MakeLogBinning(nbinsZ, zmin, zmax);
187 binLimZ=AliDielectronHelper::MakeLinBinning(nbinsZ, zmin, zmax);
190 UserProfile(histClass,name,title,valTypeP,binLimX,binLimY,binLimZ,valTypeX,valTypeY,valTypeZ,option,valTypeW);
193 //_____________________________________________________________________________
194 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
197 UInt_t valTypeX, TString option,
201 // Histogram creation 1D case with arbitraty binning
204 TVectorD *binLimX=AliDielectronHelper::MakeArbitraryBinning(binning);
205 UserProfile(histClass,name,title,valTypeP,binLimX,valTypeX,option,valTypeW);
208 //_____________________________________________________________________________
209 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
211 const TVectorD * const binsX,
212 UInt_t valTypeX/*=kNoAutoFill*/, TString option,
216 // Histogram creation 1D case with arbitraty binning X
217 // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
221 isOk&=IsHistogramOk(histClass,name);
226 if(valTypeP==kNoProfile)
227 hist=new TH1F(name,title,binsX->GetNrows()-1,binsX->GetMatrixArray());
229 TString opt=""; Double_t pmin=0., pmax=0.;
230 if(!option.IsNull()) {
231 TObjArray *arr=option.Tokenize(";");
233 opt=((TObjString*)arr->At(0))->GetString();
234 if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
235 if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
238 hist=new TProfile(name,title,binsX->GetNrows()-1,binsX->GetMatrixArray(),pmin,pmax,opt.Data());
239 // printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
242 // store variales in axes
243 UInt_t valType[20] = {0};
244 valType[0]=valTypeX; valType[1]=valTypeP;
245 StoreVariables(hist, valType);
246 hist->SetUniqueID(valTypeW); // store weighting variable
248 // store which variables are used
249 for(Int_t i=0; i<4; i++) fUsedVars->SetBitNumber(valType[i],kTRUE);
250 fUsedVars->SetBitNumber(valTypeW,kTRUE);
252 // adapt the name and title of the histogram in case they are empty
253 AdaptNameTitle(hist, histClass);
255 Bool_t isReserved=fReservedWords->Contains(histClass);
256 if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
258 UserHistogramReservedWords(histClass, hist, 999);
260 UserHistogram(histClass, hist, 999);
266 //_____________________________________________________________________________
267 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
269 const TVectorD * const binsX, const TVectorD * const binsY,
270 UInt_t valTypeX/*=kNoAutoFill*/, UInt_t valTypeY/*=0*/, TString option,
274 // Histogram creation 2D case with arbitraty binning X
275 // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
279 isOk&=IsHistogramOk(histClass,name);
285 if(valTypeP==kNoProfile) {
286 hist=new TH2F(name,title,
287 binsX->GetNrows()-1,binsX->GetMatrixArray(),
288 binsY->GetNrows()-1,binsY->GetMatrixArray());
291 TString opt=""; Double_t pmin=0., pmax=0.;
292 if(!option.IsNull()) {
293 TObjArray *arr=option.Tokenize(";");
295 opt=((TObjString*)arr->At(0))->GetString();
296 if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
297 if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
300 hist=new TProfile2D(name,title,
301 binsX->GetNrows()-1,binsX->GetMatrixArray(),
302 binsY->GetNrows()-1,binsY->GetMatrixArray());
303 ((TProfile2D*)hist)->BuildOptions(pmin,pmax,opt.Data());
304 // printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
307 // store variales in axes
308 UInt_t valType[20] = {0};
309 valType[0]=valTypeX; valType[1]=valTypeY; valType[2]=valTypeP;
310 StoreVariables(hist, valType);
311 hist->SetUniqueID(valTypeW); // store weighting variable
313 // store which variables are used
314 for(Int_t i=0; i<4; i++) fUsedVars->SetBitNumber(valType[i],kTRUE);
315 fUsedVars->SetBitNumber(valTypeW,kTRUE);
317 // adapt the name and title of the histogram in case they are empty
318 AdaptNameTitle(hist, histClass);
320 Bool_t isReserved=fReservedWords->Contains(histClass);
321 if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
323 UserHistogramReservedWords(histClass, hist, 999);
325 UserHistogram(histClass, hist, 999);
333 //_____________________________________________________________________________
334 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
336 const TVectorD * const binsX, const TVectorD * const binsY, const TVectorD * const binsZ,
337 UInt_t valTypeX/*=kNoAutoFill*/, UInt_t valTypeY/*=0*/, UInt_t valTypeZ/*=0*/, TString option,
341 // Histogram creation 3D case with arbitraty binning X
342 // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
346 isOk&=IsHistogramOk(histClass,name);
353 if(valTypeP==kNoProfile) {
354 hist=new TH3F(name,title,
355 binsX->GetNrows()-1,binsX->GetMatrixArray(),
356 binsY->GetNrows()-1,binsY->GetMatrixArray(),
357 binsZ->GetNrows()-1,binsZ->GetMatrixArray());
360 TString opt=""; Double_t pmin=0., pmax=0.;
361 if(!option.IsNull()) {
362 TObjArray *arr=option.Tokenize(";");
364 opt=((TObjString*)arr->At(0))->GetString();
365 if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
366 if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
369 hist=new TProfile3D(name,title,
370 binsX->GetNrows()-1,binsX->GetMatrixArray(),
371 binsY->GetNrows()-1,binsY->GetMatrixArray(),
372 binsZ->GetNrows()-1,binsZ->GetMatrixArray());
373 ((TProfile3D*)hist)->BuildOptions(pmin,pmax,opt.Data());
374 // printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
377 // store variales in axes
378 UInt_t valType[20] = {0};
379 valType[0]=valTypeX; valType[1]=valTypeY; valType[2]=valTypeZ; valType[3]=valTypeP;
380 StoreVariables(hist, valType);
381 hist->SetUniqueID(valTypeW); // store weighting variable
383 // store which variables are used
384 for(Int_t i=0; i<4; i++) fUsedVars->SetBitNumber(valType[i],kTRUE);
385 fUsedVars->SetBitNumber(valTypeW,kTRUE);
387 // adapt the name and title of the histogram in case they are empty
388 AdaptNameTitle(hist, histClass);
390 Bool_t isReserved=fReservedWords->Contains(histClass);
391 if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
393 UserHistogramReservedWords(histClass, hist, 999);
395 UserHistogram(histClass, hist, 999);
403 //_____________________________________________________________________________
404 void AliDielectronHistos::UserHistogram(const char* histClass, Int_t ndim, Int_t *bins, Double_t *mins, Double_t *maxs, UInt_t *vars, UInt_t valTypeW)
407 // Histogram creation 4-n dimension only with linear binning
411 isOk&=(ndim<21 && ndim>3);
412 if(!isOk) { Warning("UserHistogram","Array sizes should be between 3 and 20. Not adding Histogram to '%s'.", histClass); return; }
414 // set automatic histo name
416 for(Int_t iv=0; iv < ndim; iv++)
417 name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
418 name.Resize(name.Length()-1);
420 isOk&=IsHistogramOk(histClass,name);
424 hist=new THnD(name.Data(),"", ndim, bins, mins, maxs);
426 // store variales in axes
427 StoreVariables(hist, vars);
428 hist->SetUniqueID(valTypeW); // store weighting variable
430 // store which variables are used
431 for(Int_t i=0; i<20; i++) fUsedVars->SetBitNumber(vars[i],kTRUE);
432 fUsedVars->SetBitNumber(valTypeW,kTRUE);
434 Bool_t isReserved=fReservedWords->Contains(histClass);
436 UserHistogramReservedWords(histClass, hist, 999);
438 UserHistogram(histClass, hist, 999);
443 //_____________________________________________________________________________
444 void AliDielectronHistos::UserHistogram(const char* histClass, Int_t ndim, TObjArray *limits, UInt_t *vars, UInt_t valTypeW)
447 // Histogram creation n>3 dimension only with non-linear binning
451 isOk&=(ndim<21 && ndim>3);
452 if(!isOk) { Warning("UserHistogram","Array sizes should be between 3 and 20. Not adding Histogram to '%s'.", histClass); return; }
453 isOk&=(ndim==limits->GetEntriesFast());
456 // set automatic histo name
458 for(Int_t iv=0; iv < ndim; iv++)
459 name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
460 name.Resize(name.Length()-1);
462 isOk&=IsHistogramOk(histClass,name);
467 // get number of bins
468 for(Int_t idim=0 ;idim<ndim; idim++) {
469 TVectorD *vec = (TVectorD*) limits->At(idim);
470 bins[idim]=vec->GetNrows()-1;
473 hist=new THnD(name.Data(),"", ndim, bins, 0x0, 0x0);
476 for(Int_t idim=0 ;idim<ndim; idim++) {
477 TVectorD *vec = (TVectorD*) limits->At(idim);
478 hist->SetBinEdges(idim,vec->GetMatrixArray());
481 // store variales in axes
482 StoreVariables(hist, vars);
483 hist->SetUniqueID(valTypeW); // store weighting variable
485 // store which variables are used
486 for(Int_t i=0; i<20; i++) fUsedVars->SetBitNumber(vars[i],kTRUE);
487 fUsedVars->SetBitNumber(valTypeW,kTRUE);
489 Bool_t isReserved=fReservedWords->Contains(histClass);
491 UserHistogramReservedWords(histClass, hist, 999);
493 UserHistogram(histClass, hist, 999);
498 //_____________________________________________________________________________
499 void AliDielectronHistos::UserSparse(const char* histClass, Int_t ndim, Int_t *bins, Double_t *mins, Double_t *maxs, UInt_t *vars,
503 // THnSparse creation with linear binning
508 // set automatic histo name
510 for(Int_t iv=0; iv < ndim; iv++)
511 name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
512 name.Resize(name.Length()-1);
514 isOk&=IsHistogramOk(histClass,name);
518 hist=new THnSparseD(name.Data(),"", ndim, bins, mins, maxs);
520 // store variales in axes
521 StoreVariables(hist, vars);
522 hist->SetUniqueID(valTypeW); // store weighting variable
524 // store which variables are used
525 for(Int_t i=0; i<20; i++) fUsedVars->SetBitNumber(vars[i],kTRUE);
526 fUsedVars->SetBitNumber(valTypeW,kTRUE);
528 Bool_t isReserved=fReservedWords->Contains(histClass);
530 UserHistogramReservedWords(histClass, hist, 999);
532 UserHistogram(histClass, hist, 999);
537 //_____________________________________________________________________________
538 void AliDielectronHistos::UserSparse(const char* histClass, Int_t ndim, TObjArray *limits, UInt_t *vars, UInt_t valTypeW)
541 // THnSparse creation with non-linear binning
545 isOk&=(ndim==limits->GetEntriesFast());
548 // set automatic histo name
550 for(Int_t iv=0; iv < ndim; iv++)
551 name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
552 name.Resize(name.Length()-1);
554 isOk&=IsHistogramOk(histClass,name);
559 // get number of bins
560 for(Int_t idim=0 ;idim<ndim; idim++) {
561 TVectorD *vec = (TVectorD*) limits->At(idim);
562 bins[idim]=vec->GetNrows()-1;
565 hist=new THnSparseD(name.Data(),"", ndim, bins, 0x0, 0x0);
568 for(Int_t idim=0 ;idim<ndim; idim++) {
569 TVectorD *vec = (TVectorD*) limits->At(idim);
570 hist->SetBinEdges(idim,vec->GetMatrixArray());
573 // store variales in axes
574 StoreVariables(hist, vars);
575 hist->SetUniqueID(valTypeW); // store weighting variable
577 // store which variables are used
578 for(Int_t i=0; i<20; i++) fUsedVars->SetBitNumber(vars[i],kTRUE);
579 fUsedVars->SetBitNumber(valTypeW,kTRUE);
581 Bool_t isReserved=fReservedWords->Contains(histClass);
583 UserHistogramReservedWords(histClass, hist, 999);
585 UserHistogram(histClass, hist, 999);
590 //_____________________________________________________________________________
591 void AliDielectronHistos::UserHistogram(const char* histClass, TObject* hist, UInt_t valTypes)
594 // Add any type of user histogram
597 //special case for the calss Pair. where histograms will be created for all pair classes
598 Bool_t isReserved=fReservedWords->Contains(histClass);
600 UserHistogramReservedWords(histClass, hist, valTypes);
604 if (!IsHistogramOk(histClass,hist->GetName())) return;
605 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
606 // hist->SetDirectory(0);
608 // store variables axis
609 UInt_t valType[20] = {0};
610 // incase valTypes is given old way of extracting variables
612 valType[0]=valTypes%1000; //last three digits
613 valType[1]=valTypes/1000%1000; //second last three digits
614 valType[2]=valTypes/1000000%1000; //third last three digits
615 hist->SetUniqueID(valTypes);
618 // extract variables from axis
619 FillVarArray(hist, valType);
620 StoreVariables(hist, valType);
621 hist->SetUniqueID(valType[19]); // store weighting variable
624 classTable->Add(hist);
627 //_____________________________________________________________________________
628 void AliDielectronHistos::AddClass(const char* histClass)
631 // Add a class of histograms
632 // Several classes can be added by separating them by a ';' e.g. 'class1;class2;class3'
634 TString hists(histClass);
635 TObjArray *arr=hists.Tokenize(";");
638 while ( (o=next()) ){
639 if (fHistoList.FindObject(o->GetName())){
640 Warning("AddClass","Cannot create class '%s' it already exists.",histClass);
643 if (fReservedWords->Contains(o->GetName())){
644 Error("AddClass","Pair is a reserved word, please use another name");
647 THashList *table=new THashList;
648 table->SetOwner(kTRUE);
649 table->SetName(o->GetName());
650 fHistoList.Add(table);
655 //_____________________________________________________________________________
656 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval)
659 // Fill function 1D case
661 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
663 if (!classTable || !(hist=(TH1*)classTable->FindObject(name)) ){
664 Warning("Fill","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
670 //_____________________________________________________________________________
671 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval)
674 // Fill function 2D case
676 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
678 if (!classTable || !(hist=(TH2*)classTable->FindObject(name)) ){
679 Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
682 hist->Fill(xval,yval);
685 //_____________________________________________________________________________
686 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval, Double_t zval)
689 // Fill function 3D case
691 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
693 if (!classTable || !(hist=(TH3*)classTable->FindObject(name)) ){
694 Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
697 hist->Fill(xval,yval,zval);
700 //_____________________________________________________________________________
701 void AliDielectronHistos::FillClass(const char* histClass, Int_t nValues, const Double_t *values)
704 // Fill class 'histClass' (by name)
707 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
709 Warning("FillClass","Cannot fill class '%s' its not defined. nValues %d",histClass,nValues);
713 TIter nextHist(classTable);
715 while ( (obj=(TObject*)nextHist()) ) FillValues(obj, values);
720 //_____________________________________________________________________________
721 // void AliDielectronHistos::FillClass(const char* histClass, const TVectorD &vals)
726 // FillClass(histClass, vals.GetNrows(), vals.GetMatrixArray());
729 //_____________________________________________________________________________
730 void AliDielectronHistos::UserHistogramReservedWords(const char* histClass, const TObject *hist, UInt_t valTypes)
733 // Creation of histogram for all pair types
735 TString title(hist->GetTitle());
736 // Same Event Like Sign
737 TIter nextClass(&fHistoList);
739 while ( (l=static_cast<THashList*>(nextClass())) ){
740 TString name(l->GetName());
741 if (name.Contains(histClass)){
742 TObject *h=hist->Clone();
743 // Tobject has no function SetDirectory, didn't we need this???
744 // h->SetDirectory(0);
745 ((TH1*)h)->SetTitle(Form("%s %s",title.Data(),l->GetName()));
747 UserHistogram(l->GetName(),h,valTypes);
753 //_____________________________________________________________________________
754 void AliDielectronHistos::DumpToFile(const char* file)
757 // Dump the histogram list to a newly created root file
759 TFile f(file,"recreate");
760 fHistoList.Write(fHistoList.GetName(),TObject::kSingleKey);
764 //_____________________________________________________________________________
765 TObject* AliDielectronHistos::GetHist(const char* histClass, const char* name) const
768 // return object 'name' in 'histClass'
770 THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
771 if (!classTable) return 0x0;
772 return classTable->FindObject(name);
775 //_____________________________________________________________________________
776 TH1* AliDielectronHistos::GetHistogram(const char* histClass, const char* name) const
779 // return histogram 'name' in 'histClass'
781 return ((TH1*) GetHist(histClass, name));
784 //_____________________________________________________________________________
785 TObject* AliDielectronHistos::GetHist(const char* cutClass, const char* histClass, const char* name) const
788 // return object from list of list of histograms
789 // this function is thought for retrieving histograms if a list of AliDielectronHistos is set
792 if (!fList) return 0x0;
793 THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
795 THashList *classTable=dynamic_cast<THashList*>(h->FindObject(histClass));
796 if (!classTable) return 0x0;
797 return classTable->FindObject(name);
800 //_____________________________________________________________________________
801 TH1* AliDielectronHistos::GetHistogram(const char* cutClass, const char* histClass, const char* name) const
804 // return histogram from list of list of histograms
805 // this function is thought for retrieving histograms if a list of AliDielectronHistos is set
807 return ((TH1*) GetHist(cutClass, histClass, name));
810 //_____________________________________________________________________________
811 void AliDielectronHistos::Draw(const Option_t* option)
817 TString drawStr(option);
818 TObjArray *arr=drawStr.Tokenize(";");
823 TObjString *ostr=0x0;
827 while ( (ostr=(TObjString*)nextOpt()) ){
828 currentOpt=ostr->GetString();
829 currentOpt.Remove(TString::kBoth,'\t');
830 currentOpt.Remove(TString::kBoth,' ');
833 if ( currentOpt.Contains(testOpt.Data()) ){
834 drawClasses=currentOpt(testOpt.Length(),currentOpt.Length());
841 // Bool_t same=drawOpt.Contains("same"); //FIXME not yet implemented
846 Error("Draw","When writing to a file you have to create a canvas before opening the file!!!");
854 TIter nextClass(&fHistoList);
855 THashList *classTable=0;
856 // Bool_t first=kTRUE;
857 while ( (classTable=(THashList*)nextClass()) ){
858 //test classes option
859 if (!drawClasses.IsNull() && !drawClasses.Contains(classTable->GetName())) continue;
861 Int_t nPads = classTable->GetEntries();
862 Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
863 Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
868 canvasName.Form("c%s_%s",GetName(),classTable->GetName());
869 c=(TCanvas*)gROOT->FindObject(canvasName.Data());
870 if (!c) c=new TCanvas(canvasName.Data(),Form("%s: %s",GetName(),classTable->GetName()));
875 // if (nPads>1) gVirtualPS->NewPage();
877 if (nPads>1) c->Clear();
880 if (nCols>1||nRows>1) c->Divide(nCols,nRows);
882 //loop over histograms and draw them
883 TIter nextHist(classTable);
886 while ( (h=(TH1*)nextHist()) ){
888 if ( (h->InheritsFrom(TH2::Class())) ) drawOpt="colz";
889 if (nCols>1||nRows>1) c->cd(++iPad);
890 if ( TMath::Abs(h->GetXaxis()->GetBinWidth(1)-h->GetXaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogx();
891 if ( TMath::Abs(h->GetYaxis()->GetBinWidth(1)-h->GetYaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogy();
892 if ( TMath::Abs(h->GetZaxis()->GetBinWidth(1)-h->GetZaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogz();
893 TString histOpt=h->GetOption();
895 if (histOpt.Contains("logx")) gPad->SetLogx();
896 if (histOpt.Contains("logy")) gPad->SetLogy();
897 if (histOpt.Contains("logz")) gPad->SetLogz();
898 histOpt.ReplaceAll("logx","");
899 histOpt.ReplaceAll("logy","");
900 histOpt.ReplaceAll("logz","");
901 h->Draw(drawOpt.Data());
908 // if (gVirtualPS) delete c;
911 //_____________________________________________________________________________
912 void AliDielectronHistos::Print(const Option_t* option) const
915 // Print classes and histograms
917 TString optString(option);
919 if (optString.IsNull()) PrintStructure();
925 //_____________________________________________________________________________
926 void AliDielectronHistos::PrintStructure() const
929 // Print classes and histograms in the class to stdout
932 TIter nextClass(&fHistoList);
933 THashList *classTable=0;
934 while ( (classTable=(THashList*)nextClass()) ){
935 TIter nextHist(classTable);
937 printf("+ %s\n",classTable->GetName());
938 while ( (o=nextHist()) )
939 printf("| ->%s\n",o->GetName());
942 TIter nextCutClass(fList);
943 THashList *cutClass=0x0;
944 while ( (cutClass=(THashList*)nextCutClass()) ) {
945 printf("+ %s\n",cutClass->GetName());
946 TIter nextClass(cutClass);
947 THashList *classTable=0;
948 while ( (classTable=(THashList*)nextClass()) ){
949 TIter nextHist(classTable);
951 printf("| + %s\n",classTable->GetName());
952 while ( (o=nextHist()) )
953 printf("| | ->%s\n",o->GetName());
960 //_____________________________________________________________________________
961 void AliDielectronHistos::SetHistogramList(THashList &list, Bool_t setOwner/*=kTRUE*/)
964 // set histogram classes and histograms to this instance. It will take onwnership!
966 ResetHistogramList();
967 TString name(GetName());
968 if (name == "AliDielectronHistos") SetName(list.GetName());
971 while ( (o=next()) ){
975 list.SetOwner(kFALSE);
976 fHistoList.SetOwner(kTRUE);
978 fHistoList.SetOwner(kFALSE);
982 //_____________________________________________________________________________
983 Bool_t AliDielectronHistos::SetCutClass(const char* cutClass)
986 // Assign histogram list according to cutClass
989 if (!fList) return kFALSE;
990 ResetHistogramList();
991 THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
993 Warning("SetCutClass","cutClass '%s' not found", cutClass);
996 SetHistogramList(*h,kFALSE);
1000 //_____________________________________________________________________________
1001 Bool_t AliDielectronHistos::IsHistogramOk(const char* histClass, const char* name)
1004 // check whether the histogram class exists and the histogram itself does not exist yet
1006 Bool_t isReserved=fReservedWords->Contains(histClass);
1007 if (!fHistoList.FindObject(histClass)&&!isReserved){
1008 Warning("IsHistogramOk","Cannot create histogram. Class '%s' not defined. Please create it using AddClass before.",histClass);
1011 if (GetHist(histClass,name)){
1012 Warning("IsHistogramOk","Cannot create histogram '%s' in class '%s': It already exists!",name,histClass);
1018 // //_____________________________________________________________________________
1019 // TIterator* AliDielectronHistos::MakeIterator(Bool_t dir) const
1024 // return new TListIter(&fHistoList, dir);
1027 //_____________________________________________________________________________
1028 void AliDielectronHistos::ReadFromFile(const char* file)
1031 // Read histos from file
1034 TIter nextKey(f.GetListOfKeys());
1036 while ( (key=(TKey*)nextKey()) ){
1037 TObject *o=f.Get(key->GetName());
1038 THashList *list=dynamic_cast<THashList*>(o);
1039 if (!list) continue;
1040 SetHistogramList(*list);
1046 //_____________________________________________________________________________
1047 void AliDielectronHistos::DrawSame(const char* histName, const Option_t *opt)
1050 // Draw all histograms with the same name into one canvas
1051 // if option contains 'leg' a legend will be created with the class name as caption
1052 // if option contains 'can' a new canvas is created
1055 TString optString(opt);
1056 optString.ToLower();
1057 Bool_t optLeg=optString.Contains("leg");
1058 Bool_t optCan=optString.Contains("can");
1063 c=(TCanvas*)gROOT->FindObject(Form("c%s",histName));
1064 if (!c) c=new TCanvas(Form("c%s",histName),Form("All '%s' histograms",histName));
1069 if (optLeg) leg=new TLegend(.8,.3,.99,.9);
1072 TIter next(&fHistoList);
1073 THashList *classTable=0;
1076 while ( (classTable=(THashList*)next()) ){
1077 if ( TH1 *h=(TH1*)classTable->FindObject(histName) ){
1079 h->SetLineColor(i+1);
1080 h->SetMarkerColor(i+1);
1081 h->Draw(i>0?"same":"");
1082 if (leg) leg->AddEntry(h,classTable->GetName(),"lp");
1084 max=TMath::Max(max,h->GetMaximum());
1088 leg->SetFillColor(10);
1089 leg->SetY1(.9-i*.05);
1092 if (hFirst&&(hFirst->GetYaxis()->GetXmax()<max)){
1093 hFirst->SetMaximum(max);
1097 //_____________________________________________________________________________
1098 void AliDielectronHistos::SetReservedWords(const char* words)
1101 // set reserved words
1104 (*fReservedWords)=words;
1107 //_____________________________________________________________________________
1108 void AliDielectronHistos::StoreVariables(TObject *obj, UInt_t valType[20])
1114 if (obj->InheritsFrom(TH1::Class())) StoreVariables(static_cast<TH1*>(obj), valType);
1115 else if (obj->InheritsFrom(THnBase::Class())) StoreVariables(static_cast<THnBase*>(obj), valType);
1122 //_____________________________________________________________________________
1123 void AliDielectronHistos::StoreVariables(TH1 *obj, UInt_t valType[20])
1126 // store variables in the axis (special for TProfile3D)
1129 Int_t dim = obj->GetDimension();
1131 // dimension correction for profiles
1132 if(obj->IsA() == TProfile::Class() || obj->IsA() == TProfile2D::Class() || obj->IsA() == TProfile3D::Class()) {
1138 obj->SetUniqueID(valType[3]); // Tprofile3D variable
1140 obj->GetZaxis()->SetUniqueID(valType[2]);
1141 obj->GetZaxis()->SetName(Form("%s", AliDielectronVarManager::GetValueName(valType[2])));
1143 obj->GetYaxis()->SetUniqueID(valType[1]);
1144 obj->GetYaxis()->SetName(Form("%s", AliDielectronVarManager::GetValueName(valType[1])));
1146 obj->GetXaxis()->SetUniqueID(valType[0]);
1147 obj->GetXaxis()->SetName(Form("%s", AliDielectronVarManager::GetValueName(valType[0])));
1153 //_____________________________________________________________________________
1154 void AliDielectronHistos::StoreVariables(THnBase *obj, UInt_t valType[20])
1157 // store variables in the axis
1160 Int_t dim = obj->GetNdimensions();
1162 for(Int_t it=0; it<dim; it++) {
1163 obj->GetAxis(it)->SetUniqueID(valType[it]);
1164 obj->GetAxis(it)->SetName(Form("%s", AliDielectronVarManager::GetValueName(valType[it])));
1165 obj->GetAxis(it)->SetTitle(Form("%s %s", AliDielectronVarManager::GetValueLabel(valType[it]), AliDielectronVarManager::GetValueUnit(valType[it])));
1171 //_____________________________________________________________________________
1172 void AliDielectronHistos::FillValues(TObject *obj, const Double_t *values)
1178 if (obj->InheritsFrom(TH1::Class())) FillValues(static_cast<TH1*>(obj), values);
1179 else if (obj->InheritsFrom(THnBase::Class())) FillValues(static_cast<THnBase*>(obj), values);
1185 //_____________________________________________________________________________
1186 void AliDielectronHistos::FillValues(TH1 *obj, const Double_t *values)
1189 // fill values for TH1 inherted classes
1192 Int_t dim = obj->GetDimension();
1193 Bool_t bprf = kFALSE;
1194 // UInt_t nValues = (UInt_t) AliDielectronVarManager::kNMaxValues;
1196 UInt_t valueTypes=obj->GetUniqueID();
1197 if (valueTypes==(UInt_t)AliDielectronHistos::kNoAutoFill) return;
1198 Bool_t weight = (valueTypes!=kNoWeights);
1200 if(obj->IsA() == TProfile::Class() || obj->IsA() == TProfile2D::Class() || obj->IsA() == TProfile3D::Class())
1203 UInt_t value1=obj->GetXaxis()->GetUniqueID();
1204 UInt_t value2=obj->GetYaxis()->GetUniqueID();
1205 UInt_t value3=obj->GetZaxis()->GetUniqueID();
1206 UInt_t value4=obj->GetUniqueID(); // get profile var stored in the unique ID
1208 // ask for inclusive trigger map variables
1209 if(value1!=AliDielectronVarManager::kTriggerInclONL && value1!=AliDielectronVarManager::kTriggerInclOFF &&
1210 value2!=AliDielectronVarManager::kTriggerInclONL && value2!=AliDielectronVarManager::kTriggerInclOFF &&
1211 value3!=AliDielectronVarManager::kTriggerInclONL && value3!=AliDielectronVarManager::kTriggerInclOFF &&
1212 value4!=AliDielectronVarManager::kTriggerInclONL && value4!=AliDielectronVarManager::kTriggerInclOFF ) {
1213 // no trigger map variable selected
1216 if(!bprf && !weight) obj->Fill(values[value1]); // histograms
1217 else if(!bprf && weight) obj->Fill(values[value1], values[value4]); // weighted histograms
1218 else if(bprf && !weight) ((TProfile*)obj)->Fill(values[value1],values[value2]); // profiles
1219 else ((TProfile*)obj)->Fill(values[value1],values[value2], values[value4]); // weighted profiles
1222 if(!bprf && !weight) obj->Fill(values[value1], values[value2]); // histograms
1223 else if(!bprf && weight) ((TH2*)obj)->Fill(values[value1], values[value2], values[value4]); // weighted histograms
1224 else if(bprf && !weight) ((TProfile2D*)obj)->Fill(values[value1], values[value2], values[value3]); // profiles
1225 else ((TProfile2D*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // weighted profiles
1228 if(!bprf && !weight) ((TH3*)obj)->Fill(values[value1], values[value2], values[value3]); // histograms
1229 else if(!bprf && weight) ((TH3*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // weighted histograms
1230 else if(bprf && !weight) ((TProfile3D*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // profiles
1231 else printf(" WARNING: weighting NOT possible yet for TProfile3Ds ! \n");
1236 // fill inclusive trigger map variables
1240 for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value1],i)) obj->Fill(i); }
1243 if((value1==AliDielectronVarManager::kTriggerInclOFF && value2==AliDielectronVarManager::kTriggerInclONL) ||
1244 (value1==AliDielectronVarManager::kTriggerInclONL && value2==AliDielectronVarManager::kTriggerInclOFF) ) {
1245 for(Int_t i=0; i<30; i++) {
1246 if((UInt_t)values[value1]==BIT(i)) {
1247 for(Int_t i2=0; i2<30; i2++) {
1248 if((UInt_t)values[value2]==BIT(i2)) {
1255 else if(value1==AliDielectronVarManager::kTriggerInclONL || value1==AliDielectronVarManager::kTriggerInclOFF) {
1256 for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value1],i)) obj->Fill(i, values[value2]); }
1258 else if(value2==AliDielectronVarManager::kTriggerInclONL || value2==AliDielectronVarManager::kTriggerInclOFF) {
1259 for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value2],i)) obj->Fill(values[value1], i); }
1261 else //makes no sense
1267 } //end: trigger filling
1273 //_____________________________________________________________________________
1274 void AliDielectronHistos::FillValues(THnBase *obj, const Double_t *values)
1277 // fill values for THn inherted classes
1280 const Int_t dim = obj->GetNdimensions();
1282 UInt_t valueTypes=obj->GetUniqueID();
1283 if (valueTypes==(UInt_t)AliDielectronHistos::kNoAutoFill) return;
1284 Bool_t weight = (valueTypes!=kNoWeights);
1286 UInt_t value4=obj->GetUniqueID(); // weight variable
1289 for(Int_t it=0; it<dim; it++) fill[it] = values[obj->GetAxis(it)->GetUniqueID()];
1290 if(!weight) obj->Fill(fill);
1291 else obj->Fill(fill, values[value4]);
1297 //_____________________________________________________________________________
1298 void AliDielectronHistos::FillVarArray(TObject *obj, UInt_t *valType)
1301 // extract variables stored in the axis (special for TProfile3D)
1306 // printf(" fillvararray %s \n",obj->GetName());
1308 if (obj->InheritsFrom(TH1::Class())) {
1309 valType[0]=((TH1*)obj)->GetXaxis()->GetUniqueID();
1310 valType[1]=((TH1*)obj)->GetYaxis()->GetUniqueID();
1311 valType[2]=((TH1*)obj)->GetZaxis()->GetUniqueID();
1312 valType[3]=((TH1*)obj)->GetUniqueID(); // tprofile var stored in unique ID
1314 else if (obj->InheritsFrom(THnBase::Class())) {
1315 for(Int_t it=0; it<((THn*)obj)->GetNdimensions(); it++)
1316 valType[it]=((THn*)obj)->GetAxis(it)->GetUniqueID();
1318 valType[19]=obj->GetUniqueID(); //weights
1322 //_____________________________________________________________________________
1323 void AliDielectronHistos::AdaptNameTitle(TH1 *hist, const char* histClass) {
1326 // adapt name and title of the histogram
1329 Int_t dim = hist->GetDimension();
1330 TString currentName = hist->GetName();
1331 TString currentTitle = hist->GetTitle();
1334 Bool_t bname = (currentName.IsNull());
1335 Bool_t btitle = (currentTitle.IsNull());
1336 Bool_t bprf = kFALSE;
1337 if(hist->IsA() == TProfile::Class() || hist->IsA() == TProfile2D::Class() || hist->IsA() == TProfile3D::Class())
1341 Double_t pmin=0., pmax=0.;
1342 TString option = "", calcrange="";
1343 Bool_t bStdOpt=kTRUE;
1347 option = ((TProfile3D*)hist)->GetErrorOption();
1348 pmin = ((TProfile3D*)hist)->GetTmin();
1349 pmax = ((TProfile3D*)hist)->GetTmax();
1352 option = ((TProfile2D*)hist)->GetErrorOption();
1353 pmin = ((TProfile2D*)hist)->GetZmin();
1354 pmax = ((TProfile2D*)hist)->GetZmax();
1357 option = ((TProfile*)hist)->GetErrorOption();
1358 pmin = ((TProfile*)hist)->GetYmin();
1359 pmax = ((TProfile*)hist)->GetYmax();
1362 if(option.Contains("s",TString::kIgnoreCase)) bStdOpt=kFALSE;
1363 if(pmin!=pmax) calcrange=Form("#cbar_{%+.*f}^{%+.*f}",GetPrecision(pmin),pmin,GetPrecision(pmax),pmax);
1366 UInt_t varx = hist->GetXaxis()->GetUniqueID();
1367 UInt_t vary = hist->GetYaxis()->GetUniqueID();
1368 UInt_t varz = hist->GetZaxis()->GetUniqueID();
1369 UInt_t varp = hist->GetUniqueID();
1370 Bool_t weight = (varp!=kNoWeights);
1372 // store titles in the axis
1376 hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1378 AliDielectronVarManager::GetValueLabel(varx),
1379 AliDielectronVarManager::GetValueUnit(varx))
1381 hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1383 AliDielectronVarManager::GetValueLabel(vary),
1384 AliDielectronVarManager::GetValueUnit(vary))
1386 hist->GetZaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varz),
1388 AliDielectronVarManager::GetValueLabel(varz),
1389 AliDielectronVarManager::GetValueUnit(varz))
1392 hist->SetNameTitle(AliDielectronVarManager::GetValueName(varp),
1393 Form("%s %s%s%s%s %s",
1395 (bStdOpt ? "#LT" : "RMS("),
1396 AliDielectronVarManager::GetValueLabel(varp),
1397 (bStdOpt ? "#GT" : ")"),
1399 AliDielectronVarManager::GetValueUnit(varp))
1403 hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1405 AliDielectronVarManager::GetValueLabel(varx),
1406 AliDielectronVarManager::GetValueUnit(varx))
1408 hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1410 AliDielectronVarManager::GetValueLabel(vary),
1411 AliDielectronVarManager::GetValueUnit(vary))
1413 hist->GetZaxis()->SetTitle(Form("#%ss",histClass));
1415 hist->GetZaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varz),
1417 (bStdOpt ? "#LT" : "RMS("),
1418 AliDielectronVarManager::GetValueLabel(varz),
1419 (bStdOpt ? "#GT" : ")"),
1421 AliDielectronVarManager::GetValueUnit(varz))
1423 if(weight) hist->GetZaxis()->SetTitle(Form("%s weighted %s",
1424 AliDielectronVarManager::GetValueLabel(varp),
1425 hist->GetZaxis()->GetTitle() )
1430 hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1432 AliDielectronVarManager::GetValueLabel(varx),
1433 AliDielectronVarManager::GetValueUnit(varx))
1435 hist->GetYaxis()->SetTitle(Form("#%ss",histClass));
1437 hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1439 (bStdOpt ? "#LT" : "RMS("),
1440 AliDielectronVarManager::GetValueLabel(vary),
1441 (bStdOpt ? "#GT" : ")"),
1443 AliDielectronVarManager::GetValueUnit(vary))
1445 if(weight) hist->GetYaxis()->SetTitle(Form("%s weighted %s",
1446 AliDielectronVarManager::GetValueLabel(varp),
1447 hist->GetYaxis()->GetTitle() )
1453 // create an unique name
1457 currentName+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
1458 currentName+=Form("%s_",AliDielectronVarManager::GetValueName(vary));
1459 currentName+=Form("%s",AliDielectronVarManager::GetValueName(varz));
1460 if(bprf) currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(varp),(bStdOpt ? "avg" : "rms"));
1461 if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1464 currentName+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
1465 currentName+=Form("%s",AliDielectronVarManager::GetValueName(vary));
1466 if(bprf) currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(varz),(bStdOpt ? "avg" : "rms"));
1467 if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1470 currentName+=Form("%s",AliDielectronVarManager::GetValueName(varx));
1471 if(bprf) currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(vary),(bStdOpt ? "avg" : "rms"));
1472 if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1475 // to differentiate btw. leg and pair histos
1476 if(!strcmp(histClass,"Pair")) currentName.Prepend("p");
1477 hist->SetName(currentName.Data());
1482 Int_t AliDielectronHistos::GetPrecision(Double_t value) {
1485 // computes the precision of a double
1486 // usefull for axis ranges etc
1489 Bool_t bfnd = kFALSE;
1490 Int_t precision = 0;
1493 // printf(" value %f precision %d bfnd %d \n",TMath::Abs(value*TMath::Power(10,precision)), precision, bfnd);
1494 bfnd = ((TMath::Abs(value*TMath::Power(10,precision))/1e6 - TMath::Floor(TMath::Abs(value*TMath::Power(10,precision))/1e6)) != 0.0
1497 if(!bfnd) precision++;
1500 // printf("precision for %f found to be %d \n", value, precision);