]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/macros/MakePIDqaReport.C
o make macro compilable
[u/mrichter/AliRoot.git] / ANALYSIS / macros / MakePIDqaReport.C
1 #include "TH1.h"
2 #include "TH2.h"
3 #include "TCanvas.h"
4 #include "TFile.h"
5 #include "TLine.h"
6 #include "TList.h"
7 #include "TMath.h"
8 #include "TPaveText.h"
9 #include "TObjArray.h"
10 #include "TString.h"
11 #include "TStyle.h"
12 #include "TROOT.h"
13 #include "TSystem.h"
14 #include "TF1.h"
15 #include "TPDF.h"
16 #include "TColor.h"
17
18 #include "AliPID.h"
19
20 void SetupStyle();
21 TH2* Get2DHistogramfromList(TList *pidqalist, const char* listname, const char* histoname);
22 void AddFit(TH2* h2d);
23 void PublishCanvas(TList *qaList, const char* det, const char* name, TString nadd="");
24 void SetupPadStyle();
25
26 void LoadLibs();
27 Int_t CheckLoadLibrary(const char* library);
28
29 TCanvas *fCanvas=0x0;
30
31 /*
32
33 Example (require aliroot environment)
34
35 root.exe -l -b -q $ALICE_ROOT/ANALYSIS/macros/MakePIDqaReport.C'("PIDqa.root")'
36
37 */
38
39 void MakePIDqaReport(const char* inputFile, const char* outputFile="PIDqaReport.pdf")
40 {
41   //
42   // Make a pdf file with the efficiency report
43   //
44
45   LoadLibs();
46   SetupStyle();
47
48   TFile f(inputFile);
49   if (!f.IsOpen()){
50     printf("Could not open file '%s'\n",f.GetName());
51     return;
52   }
53
54   TList *qaList = (TList*) f.Get("PIDqa");
55   if (!qaList){
56     printf("Could not fine list 'PIDqa' in file '%s'\n",f.GetName());
57     return;
58   }
59
60   fCanvas=new TCanvas;
61
62   TPDF p(outputFile);
63
64   //
65   // Invariant mass plots
66   //
67
68
69   //
70   // Make QA info
71   //
72
73   // ITS PID
74   PublishCanvas(qaList,"ITS","hNsigmaP_ITS_%s");
75
76   // TPC PID
77   PublishCanvas(qaList,"TPC","hNsigmaP_TPC_%s");
78   //   if (man->GetCurrentPeriod()=="11h"){
79     //     PublishCanvas(qaList,"TPC","hNsigmaP_TPC_%s_Hybrid","Hybrid");
80   //     PublishCanvas(qaList,"TPC","hNsigmaP_TPC_%s_OROChigh","OROChigh");
81   //   }
82
83   // TPC PID after 3 sigma TOF cut
84   PublishCanvas(qaList,"TPC_TOF","hNsigmaP_TPC_TOF_%s");
85
86   // TOF PID
87   PublishCanvas(qaList,"TOF","hNsigmaP_TOF_%s");
88
89   // TRD PID
90   fCanvas->Divide(2,3);
91   TH2 *hLikeP_TRD_3tls_electron=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_3tls_electron");
92   TH2 *hLikeP_TRD_3tls_pion=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_3tls_pion");
93   TH2 *hLikeP_TRD_4tls_electron=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_4tls_electron");
94   TH2 *hLikeP_TRD_4tls_pion=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_4tls_pion");
95   TH2 *hLikeP_TRD_5tls_electron=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_5tls_electron");
96   TH2 *hLikeP_TRD_5tls_pion=Get2DHistogramfromList(qaList,"TRD","hLikeP_TRD_5tls_pion");
97
98   /*
99    *  cTRDnsigma[countcanvas]->cd(1);
100    *  TPaveText pt3TRD(.02,.02,.49,.52);
101    *  pt3TRD.SetTextAlign(11);
102    *  pt3TRD.SetTextSizePixels(16);
103    *  pt3TRD.AddText(Form(" TRD PID QA %s.%s.%d", first.Data(), man->GetCurrentPeriod().Data(), pass));
104    *  pt3TRD.Draw();
105    */
106   fCanvas->cd(1);
107   SetupPadStyle();
108   hLikeP_TRD_3tls_electron->Draw("colz");
109   fCanvas->cd(2);
110   SetupPadStyle();
111   hLikeP_TRD_3tls_pion->Draw("colz");
112   fCanvas->cd(3);
113   SetupPadStyle();
114   hLikeP_TRD_4tls_electron->Draw("colz");
115   fCanvas->cd(4);
116   SetupPadStyle();
117   hLikeP_TRD_4tls_pion->Draw("colz");
118   fCanvas->cd(5);
119   SetupPadStyle();
120   hLikeP_TRD_5tls_electron->Draw("colz");
121   fCanvas->cd(6);
122   SetupPadStyle();
123   hLikeP_TRD_5tls_pion->Draw("colz");
124
125   fCanvas->Update();
126   fCanvas->Clear();
127
128   // TPC Response info
129   TObjArray *qaInfo=(TObjArray*)qaList->FindObject("QAinfo");
130   TObjArray *tpcInfo=0x0;
131   if (qaInfo && (tpcInfo=(TObjArray*)qaInfo->FindObject("TPC_info"))){
132     TObjArray *tpcSplineInfo=(TObjArray*)tpcInfo->FindObject("TPC_spline_names");
133     TObjArray *tpcConfigInfo=(TObjArray*)tpcInfo->FindObject("TPC_config_info");
134     fCanvas->Divide(1,2);
135   
136     TPaveText pt(.1,.1,.9,.9,"NDC");
137     pt.SetBorderSize(1);
138     pt.SetFillColor(0);
139     pt.SetTextSizePixels(16);
140
141     if (tpcSplineInfo){
142       for (Int_t i=0; i<tpcSplineInfo->GetEntriesFast();++i) pt.AddText(tpcSplineInfo->At(i)->GetName());
143     }
144     
145     TPaveText pt2(.1,.1,.9,.9,"NDC");
146     pt2.SetBorderSize(1);
147     pt2.SetFillColor(0);
148     pt2.SetTextSizePixels(16);
149     if (tpcConfigInfo){
150       for (Int_t i=0; i<tpcConfigInfo->GetEntriesFast();++i) pt2.AddText(tpcConfigInfo->At(i)->GetName());
151     }
152
153     fCanvas->cd(1);
154     pt.Draw();
155     fCanvas->cd(2);
156     pt2.Draw();
157     
158     fCanvas->Update();
159     fCanvas->Clear();
160   }
161   
162   delete qaList;
163
164   p.Close();
165   delete fCanvas;
166 }
167
168 void SetupStyle()
169 {
170   const Int_t NCont=255;
171
172   TStyle *st = new TStyle("mystyle","mystyle");
173   gROOT->GetStyle("Plain")->Copy((*st));
174   st->SetTitleX(0.1);
175   st->SetTitleW(0.8);
176   st->SetTitleH(0.08);
177   st->SetStatX(.9);
178   st->SetStatY(.9);
179   st->SetNumberContours(NCont);
180   st->SetPalette(1,0);
181   st->SetOptStat("erm");
182   st->SetOptFit(0);
183   st->SetGridColor(kGray+1);
184   st->SetPadGridX(kTRUE);
185   st->SetPadGridY(kTRUE);
186   st->SetPadTickX(kTRUE);
187   st->SetPadTickY(kTRUE);
188   st->cd();
189
190   const Int_t NRGBs = 5;
191   Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
192   Double_t red[NRGBs]   = { 0.00, 0.00, 0.87, 1.00, 0.51 };
193   Double_t green[NRGBs] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
194   Double_t blue[NRGBs]  = { 0.51, 1.00, 0.12, 0.00, 0.00 };
195
196   TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
197
198 }
199
200 TH2* Get2DHistogramfromList(TList *pidqalist, const char* listname, const char* histoname)
201 {
202   TList *histolist = (TList *)pidqalist->FindObject(listname);
203   if (!histolist) {printf(" list not found \n");  return 0x0; }
204   TH2* histo = (TH2*)histolist->FindObject(histoname);
205   //   if (!histo) {printf(" histogram not found \n");  return 0x0; }
206   return histo;
207 }
208
209 void AddFit(TH2* h2d)
210 {
211   //
212   // Fit in slices and draw mean and sigma
213   //
214   TF1 *f1 = new TF1("f1", "gaus");
215   f1->SetRange(-1.5,1.5);
216   TObjArray aSlices;
217   h2d->FitSlicesY(f1, 0,-1, 0, "QNR",&aSlices);
218   aSlices.SetOwner(1);
219   TH1* hMean=(TH1*)aSlices.At(1);
220   TH1* hSigma=(TH1*)aSlices.At(2);
221   TH1* hChi2=(TH1*)aSlices.At(3);
222   hChi2->Scale(1./10.);
223   aSlices.AddAt(0x0,1);
224   aSlices.AddAt(0x0,2);
225   aSlices.AddAt(0x0,3);
226   hMean->SetMarkerStyle(20);
227   hMean->SetMarkerSize(0.3);
228   hMean->SetOption("same");
229   h2d->GetListOfFunctions()->Add(hMean);
230   hSigma->SetMarkerStyle(20);
231   hSigma->SetMarkerSize(0.3);
232   hSigma->SetOption("same");
233   hSigma->SetMarkerColor(kMagenta);
234   h2d->GetListOfFunctions()->Add(hSigma);
235   hChi2->SetOption("same");
236   hChi2->SetMarkerColor(kMagenta + 2);
237   hChi2->SetLineColor(kMagenta + 2);
238   h2d->GetListOfFunctions()->Add(hChi2);
239
240   TLine *l=0x0;
241   l=new TLine(h2d->GetXaxis()->GetXmin(),0,h2d->GetXaxis()->GetXmax(),0);
242   l->SetLineStyle(2);
243   h2d->GetListOfFunctions()->Add(l);
244   l=new TLine(h2d->GetXaxis()->GetXmin(),1,h2d->GetXaxis()->GetXmax(),1);
245   l->SetLineStyle(2);
246   h2d->GetListOfFunctions()->Add(l);
247 }
248
249 void PublishCanvas(TList *qaList, const char* det, const char* name, TString nadd)
250 {
251   //
252   // draw all nSigma + signal histo
253   //
254
255
256   TObjArray arrHistos;
257
258   TPaveText pt(.1,.1,.9,.9,"NDC");
259   pt.SetBorderSize(1);
260   pt.SetFillColor(0);
261   pt.SetTextSizePixels(16);
262   pt.AddText(Form("%s PID QA",det));
263   if (!nadd.IsNull()){
264     pt.AddText(nadd.Data());
265     nadd.Prepend("_");
266   }
267   arrHistos.Add(&pt);
268
269   TH2 *hSig=Get2DHistogramfromList(qaList,det,Form("hSigP_%s",det));
270   if (hSig){
271     hSig->SetOption("colz");
272     arrHistos.Add(hSig);
273   }
274
275   for (Int_t i=0;i<AliPID::kSPECIESC;++i){
276     //     for (Int_t i=0;i<AliPID::kSPECIES;++i){
277     if (i==(Int_t)AliPID::kMuon) continue;
278     TH2 *h=Get2DHistogramfromList(qaList,det,Form(name,AliPID::ParticleName(i)));
279     if (!h) continue;
280     h->SetOption("colz");
281     AddFit(h);
282     arrHistos.Add(h);
283   }
284
285   Int_t nPads=arrHistos.GetEntriesFast();
286   Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
287   Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
288
289   
290   fCanvas->Divide(nCols,nRows);
291
292
293   for (Int_t i=0; i<nPads;++i) {
294     fCanvas->cd(i+1);
295     SetupPadStyle();
296     arrHistos.At(i)->Draw();
297   }
298
299   fCanvas->Update();
300   fCanvas->Clear();
301
302 }
303
304 void SetupPadStyle()
305 {
306   gPad->SetLogx();
307   gPad->SetLogz();
308   gPad->SetGridx();
309   gPad->SetGridy();
310   gPad->SetTickx();
311   gPad->SetTicky();
312 }
313
314 void LoadLibs()
315 {
316   CheckLoadLibrary("libCore");
317   CheckLoadLibrary("libPhysics");
318   CheckLoadLibrary("libMinuit");
319   CheckLoadLibrary("libGui");
320   CheckLoadLibrary("libXMLParser");
321   
322   CheckLoadLibrary("libGeom");
323   CheckLoadLibrary("libVMC");
324   
325   CheckLoadLibrary("libNet");
326   CheckLoadLibrary("libTree");
327   CheckLoadLibrary("libProof");
328   
329   CheckLoadLibrary("libSTEERBase");
330   CheckLoadLibrary("libESD");
331   CheckLoadLibrary("libCDB");
332   CheckLoadLibrary("libRAWDatabase");
333   CheckLoadLibrary("libRAWDatarec");
334   CheckLoadLibrary("libANALYSIS");
335   CheckLoadLibrary("libSTEER");
336   
337   CheckLoadLibrary("libSTAT");
338   
339   CheckLoadLibrary("libAOD");
340   CheckLoadLibrary("libOADB");
341   CheckLoadLibrary("libANALYSISalice");
342   CheckLoadLibrary("libCORRFW");
343   
344   
345   CheckLoadLibrary("libTPCbase");
346 }
347
348 Int_t CheckLoadLibrary(const char* library)
349 {
350   // checks if a library is already loaded, if not loads the library
351   
352   if (strlen(gSystem->GetLibraries(Form("%s.so", library), "", kFALSE)) > 0)
353     return 1;
354   
355   return gSystem->Load(library);
356 }