]> git.uio.no Git - u/mrichter/AliRoot.git/blob - OADB/AliOADBPWG2Spectra.cxx
Use static method from OADBContainer instead of AnalysisManager
[u/mrichter/AliRoot.git] / OADB / AliOADBPWG2Spectra.cxx
1 #include "AliOADBPWG2Spectra.h"
2 #include "AliOADBContainer.h"
3 #include "TNamed.h"
4 #include "TString.h"
5 #include "TH1D.h"
6 #include "TObject.h"
7 #include "TList.h"
8 #include "AliAnalysisManager.h"
9 #include "TBrowser.h"
10 #include "AliLog.h"
11
12 ClassImp(AliOADBPWG2Spectra)
13
14 const char * AliOADBPWG2Spectra::fgkDetectorNames[] = {"ITS", "ITSTPC", "TPC", "TOF", "TOFTPC", "Dummy", "Dummy"};
15 const char * AliOADBPWG2Spectra::fgkPidTypeNames[]  = {"GaussFit", "NSigma", "Bayes", "Kinks"};
16 const char * AliOADBPWG2Spectra::fgkChargeTags[]    = {"Pos", "Neg"};
17 const char * AliOADBPWG2Spectra::fgkParticleNames[] = {"Pion", "Kaon", "Proton"};
18
19
20 AliOADBPWG2Spectra::AliOADBPWG2Spectra():
21 TNamed("Dummy", "OADB Object for PWG2 Spectra" ), fHistos(0)
22 {
23   // ctor
24   
25
26
27 }
28 AliOADBPWG2Spectra::AliOADBPWG2Spectra(const char* name) :
29 TNamed(name, "OADB Object for PWG2 Spectra" ), fHistos(0) 
30
31 {
32   // ctor
33   // name is appended to all histos (e.g. "Corrected")
34
35   Init();
36
37 }
38
39 AliOADBPWG2Spectra::~AliOADBPWG2Spectra() {
40   // dtor
41   if(fHistos) delete fHistos;
42 }
43
44 void AliOADBPWG2Spectra::Init() {
45   fHistos = new TList();
46 }
47
48
49 const char * AliOADBPWG2Spectra::GetOADBPWG2SpectraFileName()  {
50   // get file name to the OADB
51   static TString filename;
52   filename.Form("%s/PWG2/SPECTRA/spectraResults.root", AliOADBContainer::GetOADBPath()); 
53   return filename.Data();
54
55 }
56 const char * AliOADBPWG2Spectra::GetHistoName(Int_t det, Int_t pidType, Int_t part, 
57                                                      Int_t charge, const char * centrTag, Int_t centrBin) {
58
59   // Returns histogram name
60   // h[Name]_[Detector(s)]_[PIDType]_[Particle]_[Pos|Neg]_[MultiplicityOrCentralityIndex]
61   // where "name" is the name of this container
62
63   
64   static TString histoName;
65   if (centrTag) {
66     if(!strcmp(centrTag,"MB")){
67       // don't put a index for MB spectra
68       histoName.Form("h%s_%s_%s_%s_%s_%s", GetName(), fgkDetectorNames[det], fgkPidTypeNames[pidType], fgkParticleNames[part], fgkChargeTags[charge], centrTag);
69     }
70     else {
71       histoName.Form("h%s_%s_%s_%s_%s_%s_%d", GetName(), fgkDetectorNames[det], fgkPidTypeNames[pidType], fgkParticleNames[part], fgkChargeTags[charge], centrTag, centrBin);
72     }
73   }
74   else 
75     histoName.Form("h%s_%s_%s_%s_%s",       GetName(), fgkDetectorNames[det], fgkPidTypeNames[pidType], fgkParticleNames[part], fgkChargeTags[charge]);
76
77   return histoName.Data();
78 }
79
80 TH1D * AliOADBPWG2Spectra::GetHisto(Int_t det, Int_t pidType, Int_t part, 
81                                     Int_t charge, const char * centrTag, Int_t centrBin){
82
83   // Get an histogram from the list
84   const char * name = GetHistoName(det,pidType,part,charge,centrTag,centrBin);
85   TH1D * h = (TH1D*) fHistos->FindObject(name);
86   return h;
87
88 }
89
90 void  AliOADBPWG2Spectra::AddHisto(TH1D * h, Int_t det, Int_t pidType, Int_t part, 
91                                     Int_t charge, const char * centrTag, Int_t centrBin) {
92   // Add a histogram to the list
93   // Rename and rebinn it if necessary
94   
95   if(!h) {
96     AliWarning("Empty pointer to histogram");
97     return;
98   }
99   
100   static TH1D * htest = BookHisto(kDetDummy, kGaussFit,kPion, kPos);
101   const char * name = GetHistoName(det,pidType,part,charge,centrTag,centrBin);
102   if(!CompareBinning(h,htest)){
103     AliWarning("Histo have different binning! Rebinning to standard"){
104       h = GetHistoStandardBinning(h,det,pidType,part,charge,centrTag,centrBin);
105     }
106   }
107   if(!fHistos) {
108     AliError("fHistos not allocated!!");
109     return;
110   }
111     
112   TH1D * hold = (TH1D*) fHistos->FindObject(name);
113   if (hold) fHistos->Remove(hold);
114   delete hold;
115   if(strcmp(h->GetName(),name)){
116     AliWarning(Form("Histogram names are not consinstent %s-%s, resetting", h->GetName(),name));
117     h->SetName(name); 
118   }
119   fHistos->Add(h);
120
121
122 }
123
124 TH1D * AliOADBPWG2Spectra::BookHisto(Int_t det, Int_t pidType, Int_t part, 
125                                      Int_t charge, const char * centrTag, Int_t centrBin) {
126
127   // book a histogram according to the template. All the histograms
128   // should have the same binning (if possible/reasonable) to
129   // facilitate the compoarison and the combination
130
131   const Float_t templBins[] = {0.05,0.1,0.12,0.14,0.16,0.18,0.20,0.25,0.30,0.35,0.4,0.45,0.5,0.55,0.6,0.65,0.7,0.75,0.8,0.85,0.9,0.95,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0,2.1,2.2,2.3,2.4,2.5,2.6,2.7,2.8,2.9,3.0,3.2,3.4,3.6,3.8,4.0,4.2,4.4,4.6,4.8,5.0};
132   Int_t nbinsTempl=52;
133   const char * name = GetHistoName(det,pidType,part,charge,centrTag,centrBin);
134   TH1D * h = new TH1D(name, name, nbinsTempl, templBins);
135   return h;
136
137 }
138
139 void AliOADBPWG2Spectra::Browse(TBrowser *b)
140 {
141   // Browse this object.
142    // If b=0, there is no Browse call TObject::Browse(0) instead.
143    //         This means TObject::Inspect() will be invoked indirectly
144
145
146   if (b) {
147     b->Add(fHistos);        
148   }     
149    else
150       TObject::Browse(b);
151 }
152
153 TH1D * AliOADBPWG2Spectra::GetHistoStandardBinning(const TH1D* h, Int_t det, Int_t pidType, Int_t part, 
154                                                    Int_t charge, const char * centrTag, Int_t centrBin) {
155   // Returns a histo with the standard binning and the same content of h
156   // if the bins of h are not a subset of the standard binning, it crashes with a fatal error
157   // under and overflows are ignored
158   
159   // 1. Create a histogram with the standard binning
160   TH1D * hStd = BookHisto(det,  pidType,  part, charge, centrTag, centrBin);
161   Int_t nBinsH1=hStd->GetNbinsX();
162   Int_t nBinsH2=h->GetNbinsX();
163   // Loop over standard bins, 
164   for(Int_t iBin=1; iBin<=nBinsH1; iBin++){
165     Float_t lowPtH1 =hStd->GetBinLowEdge(iBin);
166     Float_t binWidH1=hStd->GetBinWidth(iBin);
167     // Loop over H2 bins and find overlapping bins to H1
168     for(Int_t jBin=1; jBin<=nBinsH2; jBin++){
169       Float_t lowPtH2=h->GetBinLowEdge(jBin);
170       Float_t binWidH2=h->GetBinWidth(jBin);
171       if(TMath::Abs(lowPtH1-lowPtH2)<0.001 && TMath::Abs(binWidH2-binWidH1)<0.001){
172         hStd->SetBinContent(iBin, h->GetBinContent(jBin));
173         hStd->SetBinError  (iBin, h->GetBinError  (jBin));
174         break;
175       }
176       if(TMath::Abs(lowPtH1-lowPtH2)<0.001){
177         AliError(Form("Found partially overlapping bins! [(%f,%f)(%f,%f)]",lowPtH1,binWidH1,lowPtH2,binWidH2));
178         continue;
179       }
180     }
181   }
182   return hStd;
183 }
184
185 Bool_t AliOADBPWG2Spectra::CompareBinning(TH1 * h1, TH1 * h2) {
186
187   // returns true if h1 and h2 have the same binning
188   Int_t nbins1 = h1->GetNbinsX();
189   Int_t nbins2 = h2->GetNbinsX();
190   
191   if(nbins1 != nbins2) return kFALSE;
192   
193   for(Int_t ibin = 1; ibin <= nbins1; ibin++){
194     if(TMath::Abs(h1->GetBinLowEdge(ibin) - h2->GetBinLowEdge(ibin))>0.001) return kFALSE;
195     if(TMath::Abs(h1->GetBinWidth(ibin) - h2->GetBinWidth(ibin))>0.001) return kFALSE;
196   }
197   
198   return kTRUE;
199 }