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