]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/macros/PIDCalib/compareSplines.C
- Added classes and macros for TPC PID calibration
[u/mrichter/AliRoot.git] / PWGPP / TPC / macros / PIDCalib / compareSplines.C
1 #include "TH2D.h"
2 #include "THnSparse.h"
3 #include "TCanvas.h"
4 #include "TList.h"
5 #include "TString.h"
6 #include "TFile.h"
7 #include "TGraph.h"
8 #include "TMath.h"
9 #include "TObjArray.h"
10 #include "THnSparse.h"
11 #include "TSpline.h"
12
13 #include <iostream>
14 #include <iomanip>
15
16 //__________________________________________________________________________________________
17 TSpline3* loadSplines(TFile* file, TObjArray* arr, TString name, Bool_t useOADB)
18 {
19   if (useOADB) {
20     if (!arr)
21       return 0x0;
22     
23     return (TSpline3*)arr->FindObject(name.Data());
24   }
25   
26   if (!file)
27     return 0x0;
28   
29   return (TSpline3*)file->Get(name.Data());
30 }
31
32
33
34 //__________________________________________________________________________________________
35 Int_t compareSplines(TString pathNameSplines1, TString pathNameSplines2,
36                      TString period = "LHC10D_PASS2", TString dataType1 = "DATA", TString beamType1 = "PP",
37                      TString period2 = "", TString dataType2 = "", TString beamType2 = "",
38                      Bool_t useOADBforFirstSplines = kFALSE, Bool_t useOADBforSecondSplines = kFALSE,
39                      TString displayNameSplines1 = "", TString displayNameSplines2 = "")
40
41   if (period2.IsNull())
42     period2 = period;
43   
44   if (dataType2.IsNull())
45     dataType2 = dataType1;
46   
47   if (beamType2.IsNull())
48     beamType2 = beamType1;
49   
50   if (displayNameSplines1.IsNull())
51     displayNameSplines1 = "Spline 1";
52   
53   if (displayNameSplines2.IsNull())
54     displayNameSplines2 = "Spline 2";
55   
56   TFile* f1 = 0x0;
57   f1 = TFile::Open(pathNameSplines1.Data());
58   if (!f1)  {
59     std::cout << "Failed to open file \"" << pathNameSplines1.Data() << "\"!" << std::endl;
60     return -1;
61   }
62   
63   TFile* f2 = 0x0;
64   f2 = TFile::Open(pathNameSplines2.Data());
65   if (!f2)  {
66     std::cout << "Failed to open file \"" << pathNameSplines2.Data() << "\"!" << std::endl;
67     return -1;
68   }
69
70   TCanvas* c = new TCanvas("c", "",  100,10,1380,800);
71   c->SetLogx(kTRUE);
72   
73   TH2D* hDummy = new TH2D("hDummy", Form("; p_{TPC} (GeV/c); %s / %s", displayNameSplines1.Data(), displayNameSplines2.Data()),
74                           1000, 0.15, 60, 1000, 0.9, 1.1);
75   hDummy->GetYaxis()->SetLabelSize(0.03);
76   hDummy->GetYaxis()->SetTitleSize(0.05);
77   hDummy->GetYaxis()->SetTitleOffset(1.);
78   hDummy->GetXaxis()->SetNoExponent(kTRUE);
79   hDummy->GetXaxis()->SetMoreLogLabels(kTRUE);
80   hDummy->GetXaxis()->SetTitleSize(0.05);
81   hDummy->GetXaxis()->SetLabelSize(0.04);
82   hDummy->GetXaxis()->SetTitleOffset(1.1);
83   hDummy->SetStats(kFALSE);
84   hDummy->Draw("colz");
85   
86   const Int_t nPoints = 20000;
87   const Float_t stepSize = 0.02;
88   const Float_t stepSizeEl = stepSize * 100.;
89   const Int_t nPointsEl = nPoints * 10.;
90   
91   
92   
93   TSpline3* splPion2 = 0x0;
94   TSpline3* splKaon2 = 0x0;
95   TSpline3* splElectron2 = 0x0;
96   TSpline3* splProton2 = 0x0;
97   
98   
99   TObjArray* arr = 0x0;
100   TObjArray* arr2 = 0x0;
101   
102   if (useOADBforFirstSplines) {
103     arr = (TObjArray*)f1->Get("TPCPIDResponse");
104     
105     if (!arr)  {
106       std::cout << "Failed to load first array \"TPCPIDResponse\"!" << std::endl;
107       return -1;
108     }
109   }
110   
111   if (useOADBforSecondSplines) {
112     arr2 = (TObjArray*)f2->Get("TPCPIDResponse");
113     
114     if (!arr2)  {
115       std::cout << "Failed to second array \"TPCPIDResponse\"!" << std::endl;
116       return -1;
117     }
118   }
119   
120   TSpline3* splPion = loadSplines(f1, arr, Form("TSPLINE3_%s_PION_%s_%s_MEAN", dataType1.Data(), period.Data(), beamType1.Data()), useOADBforFirstSplines);
121   splPion2 = loadSplines(f2, arr2, Form("TSPLINE3_%s_PION_%s_%s_MEAN", dataType2.Data(), period2.Data(), beamType2.Data()), useOADBforSecondSplines);
122   TGraph* gPion = new TGraph(nPoints);
123   gPion->SetTitle("#pi");
124   gPion->SetFillStyle(0);
125   gPion->SetFillColor(kWhite);
126   for (Int_t i = 0; i < nPoints; i++) {
127     gPion->SetPoint(i, (0. + i * stepSize) * 0.1396, (splPion->Eval((0. + i * stepSize)) / splPion2->Eval((0. + i * stepSize))));
128   }
129   gPion->SetLineColor(kRed);
130   gPion->SetMarkerColor(kRed);
131   gPion->Draw("same");
132
133   TSpline3* splProton = loadSplines(f1, arr, Form("TSPLINE3_%s_PROTON_%s_%s_MEAN", dataType1.Data(), period.Data(), beamType1.Data()), 
134                                     useOADBforFirstSplines);
135   splProton2 = loadSplines(f2, arr2, Form("TSPLINE3_%s_PROTON_%s_%s_MEAN", dataType2.Data(), period2.Data(), beamType2.Data()), useOADBforSecondSplines);
136   TGraph* gProton = new TGraph(nPoints);
137   gProton->SetTitle("p");
138   gProton->SetFillStyle(0);
139   gProton->SetFillColor(kWhite);
140   for (Int_t i = 0; i < nPoints; i++) {
141     gProton->SetPoint(i, (0. + i * stepSize) * 0.938, (splProton->Eval((0. + i * stepSize)) / splProton2->Eval((0. + i * stepSize))));
142   }
143   gProton->SetLineColor(kBlue);
144   gProton->SetMarkerColor(kBlue);
145   gProton->Draw("same");
146
147   TSpline3* splKaon = loadSplines(f1, arr, Form("TSPLINE3_%s_KAON_%s_%s_MEAN", dataType1.Data(), period.Data(), beamType1.Data()), useOADBforFirstSplines);
148   splKaon2 = loadSplines(f2, arr2, Form("TSPLINE3_%s_KAON_%s_%s_MEAN", dataType2.Data(), period2.Data(), beamType2.Data()), useOADBforSecondSplines);
149   TGraph* gKaon = new TGraph(nPoints);
150   gKaon->SetTitle("K");
151   gKaon->SetFillStyle(0);
152   gKaon->SetFillColor(kWhite);
153   for (Int_t i = 0; i < nPoints; i++) {
154     gKaon->SetPoint(i, (0. + i * stepSize) * 0.49368, (splKaon->Eval((0. + i * stepSize)) / splKaon2->Eval((0. + i * stepSize))));
155   }
156   gKaon->SetLineColor(kGreen);
157   gKaon->SetMarkerColor(kGreen);
158   gKaon->Draw("same");
159
160   TSpline3* splElectron = loadSplines(f1, arr, Form("TSPLINE3_%s_ELECTRON_%s_%s_MEAN", dataType1.Data(), period.Data(), beamType1.Data()), 
161                                       useOADBforFirstSplines);
162   splElectron2 = loadSplines(f2, arr2, Form("TSPLINE3_%s_ELECTRON_%s_%s_MEAN", dataType2.Data(), period2.Data(), beamType2.Data()),
163                              useOADBforSecondSplines);
164   TGraph* gElectron = new TGraph(nPoints);
165   gElectron->SetTitle("e");
166   gElectron->SetFillStyle(0);
167   gElectron->SetFillColor(kWhite);
168   for (Int_t i = 0; i < nPointsEl; i++) {
169     gElectron->SetPoint(i, (0. + i * stepSizeEl) * 0.000511, 
170                         (splElectron->Eval((0. + i * stepSizeEl)) / splElectron2->Eval((0. + i * stepSizeEl))));
171   }
172   gElectron->SetLineColor(kMagenta);
173   gElectron->SetMarkerColor(kMagenta);
174   gElectron->Draw("same");
175   
176   TLegend* leg = c->BuildLegend();
177   leg->GetListOfPrimitives()->RemoveAt(0);
178   leg->SetFillColor(kWhite);
179   c->SetGridx(kTRUE);
180   c->SetGridy(kTRUE);
181   
182   c->Update();
183   
184   return 0;
185 }