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