]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/PlotSPDpwgppQA.C
Coverity fix
[u/mrichter/AliRoot.git] / ITS / PlotSPDpwgppQA.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TSystem.h>
4 #include <Riostream.h>
5 #include <TFile.h>
6 #include <TH1F.h>
7 #include <TH2D.h>
8 #include <TH2F.h>
9 #include <TF1.h>
10 #include <TStyle.h>
11 #include <TCanvas.h>
12 #include <TMath.h>
13 #include <TROOT.h>
14 #include <TLine.h>
15 #include <TString.h>
16 #include <TPaveText.h>
17 #endif
18
19 void ratiomodules();
20 void ratiochips();
21 void mapsinner(Bool_t isShowMaps=kFALSE);
22 void phiTracklet();
23 void phiTrackletsZ();
24 void foEfficiency(); 
25
26
27 TFile *fData=0x0;
28 TFile *fMc=0x0;
29
30 TList *fListData = 0x0;
31 TList  *fListMc  = 0x0;
32
33 TString fTitleData = "";
34 TString fTitleMc = "";
35
36
37 void PlotSPDpwgppQA(TString data, TString mc, TString titleData = "[Data]", TString titleMc = "[MC]", Bool_t isGeneralTrain = kFALSE){
38
39   fTitleData=titleData;
40   fTitleMc=titleMc;
41
42   gROOT->SetStyle("Plain");
43   gStyle->SetPalette(1);
44   gStyle->SetOptStat(0);
45   gStyle->SetOptFit(111);
46
47  fData = TFile::Open(data.Data());
48   if(!isGeneralTrain) fListData = (TList*)fData->Get("chist");
49   else {
50     TDirectoryFile *spddata = (TDirectoryFile*)fData->Get("SPD_Performance");
51     spddata->cd();
52     fListData = (TList*)spddata->Get("coutput1");
53   }
54
55   Double_t nevtsData = ((TH1I*)(fListData->FindObject("hEventsProcessed")))->GetEntries();
56   printf("   #events in %s : %f \n",fTitleData.Data(),nevtsData);
57
58   fMc = TFile::Open(mc.Data());
59   if(!isGeneralTrain) fListMc = (TList*)fMc->Get("chist");
60   else {
61     TDirectoryFile *spdmc = (TDirectoryFile*)fMc->Get("SPD_Performance");
62     spdmc->cd();
63     fListMc = (TList*)spdmc->Get("coutput1");
64   }
65   Double_t nevtsMc = ((TH1I*)(fListMc->FindObject("hEventsProcessed")))->GetEntries();
66   // phi projection
67
68   printf("   #events in %s : %f \n",fTitleMc.Data(),nevtsMc);
69   printf("Available functions : \n - ratiomodules() \n - ratiochips() \n - mapsinner(isShowMaps) \n - phiTracklet() \n - phiTrackletsZ() \n - foEfficiency() \n");
70   phiTracklet();
71 }
72
73 void ratiomodules(){
74
75   TH1F *data_module = ((TH1F*)fListData->FindObject("hClusterModYield"));
76   data_module->Sumw2();
77
78   TH1F *mc_module = ((TH1F*)fListMc->FindObject("hClusterModYield"));
79   mc_module->Sumw2();
80
81   TCanvas *c1 = new TCanvas("c1","c1",1200,600);
82   c1->Divide(2,1);
83   c1->cd(1);
84   data_module->SetTitle(Form("cluster yield in %s",fTitleData.Data()));
85   data_module->DrawCopy();
86   c1->cd(2);
87   mc_module->SetTitle(Form("cluster yield in %s",fTitleMc.Data()));
88   mc_module->DrawCopy();
89
90
91   TH1D *ratiomodule  = new TH1D("ratiomod",Form("Module cluster ratio %s / %s  - scaling : Integral",fTitleData.Data(),fTitleMc.Data()),mc_module->GetNbinsX(),mc_module->GetXaxis()->GetXmin(),mc_module->GetXaxis()->GetXmax());
92   ratiomodule->GetXaxis()->SetTitle("module number");
93   printf("data_module %f     -  mc_module %f    \n",data_module->GetEntries(),mc_module->GetEntries());
94   ratiomodule->Divide(data_module,mc_module,mc_module->Integral(),data_module->Integral(),"e");
95
96   TCanvas *ratiofit = new TCanvas("ratiofit","ratiofit",1200,600);
97   ratiofit->cd();
98  
99   TF1 *f1inner = new TF1("f1inner", "pol0", 0, 79);
100   f1inner->SetLineColor(kBlue);
101   TF1 *f1outer = new TF1("f1outer", "pol0", 80, 223);
102   f1outer->SetLineColor(kGreen);
103   
104   printf("------------- Fitting standard normalized ratio  ----------------\n");
105   ratiomodule->Fit("f1inner", "R+"); 
106   ratiomodule->Fit("f1outer", "R+");
107   ratiomodule->GetListOfFunctions()->Print();
108    
109   f1inner->Draw("same");
110   f1outer->Draw("same");
111   
112   TPaveText *parsFitEntries = new TPaveText(0.3,0.23,0.6,0.55,"NDC"); 
113   parsFitEntries->AddText("Inner ");
114   parsFitEntries->AddText(Form("#chi2 / ndf    %4.3f / %2.1i",f1inner->GetChisquare(),f1inner->GetNDF()));
115   parsFitEntries->AddText(Form("Fit :    %4.4f +- %4.4f",f1inner->GetParameter(0),f1inner->GetParError(0)));
116   parsFitEntries->AddText("");
117   parsFitEntries->AddText("Outer ");
118   parsFitEntries->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",f1outer->GetChisquare(),f1outer->GetNDF()));
119   parsFitEntries->AddText(Form("Fit :    %4.4f +- %4.4f",f1outer->GetParameter(0),f1outer->GetParError(0)));
120   ratiomodule->GetListOfFunctions()->Add(parsFitEntries);
121
122
123
124
125   TCanvas *ratiofitEvents = new TCanvas("ratiofitEvents","ratiofitEvents",1200,600);
126
127   ratiofitEvents->cd();
128   TH1D *ratiomoduleEvents  = new TH1D("ratiomodEvents",Form("Module cluster ratio %s / %s  - scaling : # events",fTitleData.Data(),fTitleMc.Data()),mc_module->GetNbinsX(),mc_module->GetXaxis()->GetXmin(),mc_module->GetXaxis()->GetXmax());
129   ratiomoduleEvents->GetXaxis()->SetTitle("module number");
130   Double_t nEvData = ((TH1I*)(fListData->FindObject("hEventsProcessed")))->GetEntries();
131   Double_t nEvMc= ((TH1I*)(fListMc->FindObject("hEventsProcessed")))->GetEntries();
132   ratiomoduleEvents->Divide(data_module,mc_module,nEvMc,nEvData,"e");
133   ratiomoduleEvents->Draw();
134
135   TF1 *fInner = new TF1("fInner", "pol0", 0, 79);
136   fInner->SetLineColor(kBlue);
137   TF1 *fOuter = new TF1("fOuter", "pol0", 80, 223);
138   fOuter->SetLineColor(kGreen);
139   printf("------------- Fitting #evts normalized ratio  ----------------\n");
140   ratiomoduleEvents->Fit("fInner", "R");
141   ratiomoduleEvents->Fit("fOuter", "R+");
142   fInner->Draw("same");
143   fOuter->Draw("same");
144   
145   
146   
147   TPaveText *parsFitEvents = new TPaveText(0.3,0.23,0.6,0.55,"NDC"); 
148   parsFitEvents->AddText("Inner ");
149   parsFitEvents->AddText(Form("#chi2 / ndf    %4.3f / %2.1i",fInner->GetChisquare(),fInner->GetNDF()));
150   parsFitEvents->AddText(Form("Fit :    %4.4f +- %4.4f",fInner->GetParameter(0),fInner->GetParError(0)));
151   parsFitEvents->AddText("");
152   parsFitEvents->AddText("Outer ");
153   parsFitEvents->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",fOuter->GetChisquare(),fOuter->GetNDF()));
154   parsFitEvents->AddText(Form("Fit :    %4.4f +- %4.4f",fOuter->GetParameter(0),fOuter->GetParError(0)));
155   ratiomoduleEvents->GetListOfFunctions()->Add(parsFitEvents);
156   
157
158 }
159 void mapsinner(Bool_t isShowMaps){
160   TH2D *data_mapL1 = (TH2D*)fListData->FindObject("hLocalMapL1");
161   // ------------ phi projection ---------------------
162   TH1D *data_projyL1 = data_mapL1->ProjectionY();
163
164   TString titleDataL1  = data_projyL1->GetTitle();
165   titleDataL1+=fTitleData.Data();
166   data_projyL1->SetTitle(titleDataL1.Data());
167   //data_projyL1->Rebin(10);
168   data_projyL1->SetYTitle(Form("entries / %1.2f cm",data_projyL1->GetBinWidth(0)));
169
170   TH2D *data_mapL2 = (TH2D*)fListData->FindObject("hLocalMapL2");
171   TH1D *data_projyL2 = data_mapL2->ProjectionY();
172   TString titleDataL2  = data_projyL2->GetTitle();
173   titleDataL2+=fTitleData.Data();
174   data_projyL2->SetTitle(titleDataL2.Data());
175   //data_projyL2 ->Rebin(10);
176   data_projyL2->SetYTitle(Form("entries / %1.2f cm",data_projyL2->GetBinWidth(0)));
177
178   // ------- z projection Data -----------
179   TH1D *data_projyL1z = data_mapL1->ProjectionX();
180   TString titleDataL1z  = data_projyL1z->GetTitle();
181   titleDataL1z+=fTitleData.Data();
182   data_projyL1z->SetTitle(titleDataL1z.Data());
183   //data_projyL1z->Rebin(10);
184   data_projyL1z->SetYTitle(Form("entries / %1.2f cm",data_projyL1z->GetBinWidth(0)));
185   TH1D *data_projyL2z = data_mapL2->ProjectionX();
186   TString titleDataL2z  = data_projyL2z->GetTitle();
187   titleDataL2z+=fTitleData.Data();
188   data_projyL2z->SetTitle(titleDataL2.Data());
189   //data_projyL2z ->Rebin(10);
190   data_projyL2z->SetYTitle(Form("entries / %1.2f cm",data_projyL2z->GetBinWidth(0)));
191
192   // ------------ phi projection MC ---------------------
193   TH2D *mc_mapL1 = (TH2D*)fListMc->FindObject("hLocalMapL1");
194   TH1D *mc_projyL1 = mc_mapL1->ProjectionY();
195   TString titleMCL1  = mc_projyL1->GetTitle();
196   titleMCL1+=fTitleMc.Data();
197   mc_projyL1->SetTitle(titleMCL1.Data());
198   //mc_projyL1->Rebin(10);
199   mc_projyL1->SetYTitle(Form("entries / %1.2f cm",mc_projyL1->GetBinWidth(0)));
200
201   TH2D *mc_mapL2 = (TH2D*)fListMc->FindObject("hLocalMapL2");
202   TH1D *mc_projyL2 = mc_mapL2->ProjectionY();
203   TString titleMCL2  = mc_projyL2->GetTitle();
204   titleMCL2+=fTitleMc.Data();
205   mc_projyL2->SetTitle(titleMCL2.Data());
206   //mc_projyL2->Rebin(10);
207   mc_projyL2->SetYTitle(Form("entries / %1.2f cm",mc_projyL2->GetBinWidth(0)));
208
209   // ------- z projection MC -----------
210
211   TH1D *mc_projyL1z = mc_mapL1->ProjectionX();
212   TString titleMCL1z  = mc_projyL1z->GetTitle();
213   titleMCL1z+=fTitleMc.Data();
214   mc_projyL1z->SetTitle(titleMCL1z.Data());
215   //mc_projyL1z->Rebin(10);
216   mc_projyL1z->SetYTitle(Form("entries / %1.2f cm",mc_projyL1z->GetBinWidth(0)));
217
218   TH1D *mc_projyL2z = mc_mapL2->ProjectionX();
219   TString titleMCL2z  = mc_projyL2z->GetTitle();
220   titleMCL2z+="  [ MC ]";
221   mc_projyL2z->SetTitle(titleMCL2z.Data());
222   //mc_projyL2z->Rebin(10);
223   mc_projyL2z->SetYTitle(Form("entries / %1.2f cm",mc_projyL2z->GetBinWidth(0)));
224
225   if(isShowMaps) {
226
227
228     TCanvas *cmapData = new TCanvas("cmapData","cmapData",1200,600);
229     cmapData->Divide(2,1);
230     cmapData->cd(1);
231     TString titledata1 = data_mapL1->GetTitle();
232     titledata1+=fTitleData.Data();
233     data_mapL1->SetTitle(titledata1.Data());
234     data_mapL1->Draw("colz");
235     cmapData->cd(2);
236     TString titledata2 = data_mapL2->GetTitle();
237     titledata2+=fTitleData.Data();
238     data_mapL2->SetTitle(titledata2.Data());
239     data_mapL2->Draw("colz");
240
241     TCanvas *cmapMc = new TCanvas("cmapMc","cmapMc",1200,600);
242     cmapMc->Divide(2,1);
243
244     cmapMc->cd(1);
245     TString titlemc1 = mc_mapL1->GetTitle();
246     titlemc1+=fTitleMc.Data();
247     mc_mapL1->SetTitle(titlemc1.Data());
248     mc_mapL1->Draw("colz");
249     cmapMc->cd(2);
250     TString titlemc2 = mc_mapL2->GetTitle();
251     titlemc2+=fTitleMc.Data();
252     mc_mapL2->SetTitle(titlemc2.Data());
253     mc_mapL2->Draw("colz");
254   }
255   // booking ratios 
256   // projection phi
257
258   TH1D *ratioL1 = new TH1D("ratioL1","Data / MC - Layer 1",mc_projyL1->GetNbinsX(),mc_projyL1->GetXaxis()->GetXmin(),mc_projyL1->GetXaxis()->GetXmax());
259   ratioL1->SetXTitle("(NOT GLOBAL)       direction #varphi [cm]");
260   TH1D *ratioL2 = new TH1D("ratioL2","Data / MC - Layer 2",mc_projyL2->GetNbinsX(),mc_projyL2->GetXaxis()->GetXmin(),mc_projyL2->GetXaxis()->GetXmax());
261   ratioL2->SetXTitle("(NOT GLOBAL)       direction #varphi [cm]");
262   // projection z 
263  
264   TH1D *ratioL1z = new TH1D("ratioL1z","Data / MC - Layer 1",mc_projyL1z->GetNbinsX(),mc_projyL1z->GetXaxis()->GetXmin(),mc_projyL1z->GetXaxis()->GetXmax());
265   ratioL1z->SetXTitle("(NOT GLOBAL)       direction z [cm]");
266   TH1D *ratioL2z = new TH1D("ratioL2z","Data / MC - Layer 2",mc_projyL2z->GetNbinsX(),mc_projyL2z->GetXaxis()->GetXmin(),mc_projyL2z->GetXaxis()->GetXmax());
267   ratioL2->SetXTitle("(NOT GLOBAL)       direction z [cm]");
268
269   // making the ratios
270   ratioL1->Divide(data_projyL1,mc_projyL1,mc_projyL1->GetEntries(),data_projyL1->GetEntries(),"e");
271   ratioL2->Divide(data_projyL2,mc_projyL2,mc_projyL2->GetEntries(),data_projyL2->GetEntries(),"e");
272
273   ratioL1z->Divide(data_projyL1z,mc_projyL1z,mc_projyL1z->GetEntries(),data_projyL1z->GetEntries(),"e");
274   ratioL2z->Divide(data_projyL2z,mc_projyL2z,mc_projyL2z->GetEntries(),data_projyL2z->GetEntries(),"e");
275
276   TCanvas *c = new TCanvas("cRatioPhi","cRatioPhi",1200,600);
277   c->Divide(2,1);
278   c->cd(1);
279   ratioL1->SetYTitle(Form("ratio / %2.2f cm",ratioL1->GetBinWidth(0)));
280   ratioL1->Draw();
281   c->cd(2);
282   ratioL2->SetYTitle(Form("ratio / %2.2f cm",ratioL2->GetBinWidth(0)));
283   ratioL2->Draw();
284
285   //TLine *l1z = new TLine(-16.5,1,16.5,1);
286   //TLine *l2z = new TLine(-16.5,1,16.5,1);
287   TCanvas *cz = new TCanvas("cRatioZ","cRatioZ",1200,600);
288   cz->Divide(2,1);
289   cz->cd(1);
290   ratioL1z->SetYTitle(Form("ratio / %2.2f cm",ratioL1z->GetBinWidth(0)));
291   ratioL1z->Draw();
292   cz->cd(2);
293   ratioL2z->SetYTitle(Form("ratio / %2.2f cm",ratioL2z->GetBinWidth(0)));
294   ratioL2z->Draw();
295
296 }
297 void phiTracklet(){
298
299   TCanvas *rawDist = new TCanvas("rawDist"," raw distributions ",1200,800);
300   rawDist->Divide(2,2);
301
302   TH2F *trackData = (TH2F*)fListData->FindObject("hSPDphivsSPDeta");
303   trackData->SetTitle(Form("%s %s",trackData->GetTitle(),fTitleData.Data()));
304   TH1D *trackDataPhi = trackData->ProjectionY();
305   if(!trackDataPhi) printf("NO 1 \n");
306   //trackDataPhi->SetTitle(Form("%s %s",trackDataPhi->GetTitle(),fTitleData.Data()));
307   rawDist->cd(1);
308   trackDataPhi->SetLineColor(kRed);
309   trackDataPhi->DrawCopy();
310   TH1D *trackDataEta = trackData->ProjectionX();
311   if(!trackDataEta) printf("NO 2 \n");
312   //trackDataEta->SetTitle(Form("%s %s",trackDataEta->GetTitle(),fTitleData.Data()));
313   rawDist->cd(2);
314   trackDataEta->SetLineColor(kRed);
315   trackDataEta->DrawCopy();
316
317   TH1F etaData, phiData;
318   trackDataEta->Copy(etaData);
319   trackDataPhi->Copy(phiData);
320
321   TH1F etaFrac, phiFrac, mcEta, mcPhi;
322   trackDataEta->Copy(etaFrac);
323   trackDataPhi->Copy(phiFrac);
324
325   TH2F *trackMc = (TH2F*)fListMc->FindObject("hSPDphivsSPDeta");
326   trackMc->SetTitle(Form("%s %s",trackMc->GetTitle(),fTitleMc.Data()));
327
328   TCanvas *tracklets = new TCanvas("tracklets","tracklets",1200,600);
329   tracklets->Divide(2,1);
330   tracklets->cd(1);
331   tracklets->cd(1)->SetRightMargin(0.15);
332   //trackData->SetTitle(Form("%s %s",trackData->GetTitle(),fTitleData.Data()));
333   trackData->DrawCopy("colz");
334   tracklets->cd(2);
335   tracklets->cd(2)->SetRightMargin(0.15);
336   //trackMc->SetTitle(Form("%s %s",trackMc->GetTitle(),fTitleMc.Data()));
337   TH1D *h = (TH1D*)trackMc->DrawCopy("colz");
338   fTitleData.ReplaceAll(" ","");
339   fTitleMc.ReplaceAll(" ","");
340   tracklets->SaveAs(Form("trackletsPhiEtaMaps_%s_%s.png",fTitleData.Data(),fTitleMc.Data()));
341
342   TH1D *trackMcPhi = trackMc->ProjectionY();
343   trackMcPhi->SetTitle(Form("%s",h->GetTitle()));
344   rawDist->cd(3);
345   trackMcPhi->DrawCopy();
346   TH1D *trackMcEta = trackMc->ProjectionX();
347   trackMcEta->SetTitle(Form("%s",h->GetTitle()));
348   rawDist->cd(4);
349   trackMcEta->DrawCopy();
350
351   rawDist->SaveAs(Form("trackletsPhiEtaRaw_%s_%s.png",fTitleData.Data(),fTitleMc.Data()));
352
353   TH1F etaMc, phiMc;
354   trackMcEta->Copy(etaMc);
355   trackMcPhi->Copy(phiMc);
356
357   trackMcEta->Copy(mcEta);
358   trackMcPhi->Copy(mcPhi);
359
360   etaFrac.Scale(1./etaFrac.GetEntries());
361   mcEta.Scale(1./mcEta.GetEntries());
362   etaFrac.Add(&mcEta,-1);
363   etaFrac.Divide(&mcEta);
364
365   phiFrac.Scale(1./phiFrac.GetEntries());
366   mcPhi.Scale(1./mcPhi.GetEntries());
367   phiFrac.Add(&mcPhi,-1);
368   phiFrac.Divide(&mcPhi);
369
370
371   TCanvas *track = new TCanvas("track","track",1200,600);
372   track->Divide(2,1);
373   track->cd(1);
374   phiData.SetLineColor(kRed);
375   phiData.SetLineWidth(2);
376   phiData.Scale(1./phiData.GetEntries());
377   phiData.DrawCopy();
378   phiMc.Scale(1./phiMc.GetEntries());
379   phiMc.DrawCopy("same");
380   track->cd(2);
381   etaData.SetLineColor(kRed);
382   etaData.SetLineWidth(2);
383   etaData.Scale(1./etaData.GetEntries());
384   etaData.DrawCopy();
385   etaMc.Scale(1./etaMc.GetEntries());
386   etaMc.DrawCopy("same");
387   track->SaveAs(Form("trackletsPhiEtaNorm_%s_%s.png",fTitleData.Data(),fTitleMc.Data()));
388
389   TCanvas *frac = new TCanvas("frac","fractions",1200,600);
390   frac->Divide(2,1);
391   frac->cd(1);
392   etaFrac.SetTitle(Form(" #Delta#eta/#eta_{%s}   %s - %s ",fTitleMc.Data(),fTitleData.Data(),fTitleMc.Data()));
393   etaFrac.SetLineColor(1);
394   etaFrac.DrawCopy();
395   frac->cd(2);
396   phiFrac.SetTitle(Form(" #Delta#varphi/#varphi_{%s}   %s - %s ",fTitleMc.Data(),fTitleData.Data(),fTitleMc.Data()));
397   phiFrac.SetLineColor(1);
398   phiFrac.DrawCopy();
399
400   frac->SaveAs(Form("relativeRatios_%s_%s.png",fTitleData.Data(),fTitleMc.Data()));
401
402 }
403
404 void foEfficiency(){
405   TH2F *firedFoData = (TH2F*)fListData->FindObject("hFOgoodPerBCmod4");
406   TH2F *firedChipsData = (TH2F*)fListData->FindObject("hFiredChipsPerBCmod4");
407   TH2F mapBCmod;
408   firedFoData->Copy(mapBCmod); 
409   mapBCmod.Divide(firedChipsData);
410   mapBCmod.SetTitle(Form("FO eff per BCmod4 in %s ",fTitleData.Data()));
411   mapBCmod.GetYaxis()->SetNdivisions(4,kFALSE);
412   TCanvas *c = new TCanvas("mapFo"," FO eff map",800,800);
413   c->cd();
414   c->cd()->SetGridy();
415   mapBCmod.DrawCopy("colz");
416
417   TH1F bcmod[4];
418
419   TH1F *hbc = new TH1F("bc","bc",firedFoData->GetNbinsX(),(firedFoData->GetXaxis())->GetXmin(),(firedFoData->GetXaxis())->GetXmax());
420
421   for(Int_t bc=0; bc<4; bc++){
422     hbc->Clear();
423     hbc->Reset();
424     for(Int_t iBin=0; iBin<firedFoData->GetNbinsX(); iBin++){
425       hbc->SetBinContent(iBin+1,mapBCmod.GetBinContent(iBin+1,bc+1)); 
426     }
427     hbc->Copy(bcmod[bc]);
428     bcmod[bc].SetLineColor(bc+1); 
429   }
430
431
432   TH1F *h[4];
433   TCanvas *bceff = new TCanvas("bceff");
434   bceff->Divide(2,2);
435   for(Int_t iPad=0; iPad<4; iPad++){
436    TVirtualPad *pad = bceff->cd(iPad+1);
437     if(iPad<2){
438       pad->Divide(2,1);
439       Int_t idx = -1;
440       if(iPad==0) idx=0;
441       if(iPad==1) idx=2;  
442       pad->cd(1); bcmod[idx].DrawCopy();
443       pad->cd(2); bcmod[idx+1].DrawCopy();
444     }
445
446     if(iPad==2){
447       h[0]= (TH1F*)(bcmod[0].DrawCopy());
448       for(Int_t bb=1; bb<4; bb++){
449         h[bb] = (TH1F*)(bcmod[bb].DrawCopy("same"));
450       }
451     }
452
453     TH1F *hdiff[3];
454
455     if(iPad==3){
456       hdiff[0]=(TH1F*)(h[1]->Clone());   
457       hdiff[0]->Add(h[0],-1);
458       hdiff[0]->DrawCopy();
459       hdiff[1]=(TH1F*)(h[2]->Clone());   
460       hdiff[1]->Add(h[0],-1);
461       hdiff[1]->DrawCopy("same");
462       hdiff[2]=(TH1F*)(h[3]->Clone());   
463       hdiff[2]->Add(h[0],-1);
464       hdiff[2]->DrawCopy("same");
465     }
466   }
467 }
468
469 void phiTrackletsZ(){
470
471   TH1F* phiZposData = (TH1F*)fListData->FindObject("hSPDphiZpos");
472   phiZposData->SetLineColor(kRed);
473   phiZposData->SetLineWidth(2);
474   TH1F* phiZnegData = (TH1F*)fListData->FindObject("hSPDphiZneg");
475   phiZnegData->SetLineColor(kRed);
476   phiZnegData->SetLineWidth(2);
477   TCanvas *cZ = new TCanvas("cZ","cZ",1000,600);
478
479   cZ->Divide(2,1);
480   cZ->cd(1);
481   phiZposData->Scale(1./phiZposData->GetEntries());
482   phiZposData->DrawCopy();
483   cZ->cd(2);
484   phiZnegData->Scale(1./phiZnegData->GetEntries());
485   phiZnegData->DrawCopy();
486
487   TH1F* phiZposMc = (TH1F*)fListMc->FindObject("hSPDphiZpos");
488   TH1F* phiZnegMc = (TH1F*)fListMc->FindObject("hSPDphiZneg");
489   cZ->cd(1);
490   phiZposMc->Scale(1./phiZposMc->GetEntries());
491   phiZposMc->DrawCopy("same");
492   cZ->cd(2);
493   phiZnegMc->Scale(1./phiZnegMc->GetEntries());
494   phiZnegMc->DrawCopy("same");
495
496
497
498 void ratiochips(){
499
500   TH1F *data_chip = ((TH1F*)fListData->FindObject("hClusterYield"));
501   data_chip->Sumw2();
502
503   TH1F *mc_chip = ((TH1F*)fListMc->FindObject("hClusterYield"));
504   mc_chip->Sumw2();
505
506   TCanvas *c1chip = new TCanvas("c1chip","c1chip",1200,600);
507   c1chip->Divide(2,1);
508   c1chip->cd(1);
509   data_chip->SetTitle(Form("chip cluster yield in %s",fTitleData.Data()));
510   data_chip->DrawCopy();
511   c1chip->cd(2);
512   mc_chip->SetTitle(Form("chip cluster yield in %s",fTitleMc.Data()));
513   mc_chip->DrawCopy();
514
515
516   TH1D *ratiochip  = new TH1D("ratiochip",Form("Chip cluster ratio %s / %s  - scaling : Integral",fTitleData.Data(),fTitleMc.Data()),mc_chip->GetNbinsX(),mc_chip->GetXaxis()->GetXmin(),mc_chip->GetXaxis()->GetXmax());
517   ratiochip->GetXaxis()->SetTitle("chip number");
518   printf("data_chip %f     -  mc_chip %f    \n",data_chip->GetEntries(),mc_chip->GetEntries());
519   ratiochip->Divide(data_chip,mc_chip,mc_chip->Integral(),data_chip->Integral(),"e");
520   ratiochip->SetMarkerStyle(20);
521   ratiochip->SetMarkerSize(0.7);
522   ratiochip->SetMarkerColor(2);
523
524   TCanvas *ratiofitchip = new TCanvas("ratiofitchip","ratiofit",1200,600);
525   ratiofitchip->cd();
526   TF1 *f1innerChip = new TF1("f1innerChip", "pol0", 0, 399);
527   f1innerChip->SetLineColor(kBlue);
528   TF1 *f1outerChip = new TF1("f1outerChip", "pol0", 400, 1199);
529   f1outerChip->SetLineColor(kGreen);
530   
531   printf("------------- Fitting standard normalized ratio  ----------------\n");
532   ratiochip->Fit("f1innerChip", "R");
533   ratiochip->Fit("f1outerChip", "R+");
534   //f1inner->GetParameters(&par[0]);
535   //f1outer->GetParameters(&par[1]);
536   //f1->SetParameters(par);
537   //ratiomodule->Fit(f1, "R+");
538   f1innerChip->Draw("same");
539   f1outerChip->Draw("same");
540
541   TPaveText *parsFitChip= new TPaveText(0.3,0.23,0.6,0.55,"NDC"); 
542   parsFitChip->AddText("Inner ");
543   parsFitChip->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",f1innerChip->GetChisquare(),f1innerChip->GetNDF()));
544   parsFitChip->AddText(Form("Fit :    %4.4f +- %4.4f",f1innerChip->GetParameter(0),f1innerChip->GetParError(0)));
545   parsFitChip->AddText("");
546   parsFitChip->AddText("Outer ");
547   parsFitChip->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",f1outerChip->GetChisquare(),f1outerChip->GetNDF()));
548   parsFitChip->AddText(Form("Fit :    %4.4f +- %4.4f",f1outerChip->GetParameter(0),f1outerChip->GetParError(0)));
549   ratiochip->GetListOfFunctions()->Add(parsFitChip);
550
551
552
553   TCanvas *ratiofitChipEvents = new TCanvas("ratiofitChipEvents","ratiofitEvents",1200,600);
554
555   ratiofitChipEvents->cd();
556   TH1D *ratiochipEvents  = new TH1D("ratiochipEvents",Form("Chip cluster ratio %s / %s  - scaling : # events",fTitleData.Data(),fTitleMc.Data()),mc_chip->GetNbinsX(),mc_chip->GetXaxis()->GetXmin(),mc_chip->GetXaxis()->GetXmax());
557   ratiochipEvents->GetXaxis()->SetTitle("chip number");
558   Double_t nEvData = ((TH1I*)(fListData->FindObject("hEventsProcessed")))->GetEntries();
559   Double_t nEvMc= ((TH1I*)(fListMc->FindObject("hEventsProcessed")))->GetEntries();
560   ratiochipEvents->Divide(data_chip,mc_chip,nEvMc,nEvData,"e");
561   ratiochipEvents->Draw();
562   ratiochipEvents->SetMarkerStyle(20);
563   ratiochipEvents->SetMarkerSize(0.7);
564   ratiochipEvents->SetMarkerColor(2);
565
566   TF1 *fInnerChip = new TF1("fInnerChip", "pol0", 0, 399);
567   fInnerChip->SetLineColor(kBlue);
568   TF1 *fOuterChip = new TF1("fOuterChip", "pol0", 400, 1199);
569   fOuterChip->SetLineColor(kGreen);
570   ratiochipEvents->Fit("fInnerChip", "R");
571   ratiochipEvents->Fit("fOuterChip", "R+");
572   fInnerChip->Draw("same");
573   fOuterChip->Draw("same");
574   
575   
576   TPaveText *parsFit= new TPaveText(0.3,0.23,0.6,0.55,"NDC"); 
577   parsFit->AddText("Inner ");
578   parsFit->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",fInnerChip->GetChisquare(),fInnerChip->GetNDF()));
579   parsFit->AddText(Form("Fit :    %4.4f +- %4.4f",fInnerChip->GetParameter(0),fInnerChip->GetParError(0)));
580   parsFit->AddText("");
581   parsFit->AddText("Outer ");
582   parsFit->AddText(Form("#chi2 / ndf    %4.3f / %3.1i",fOuterChip->GetChisquare(),fOuterChip->GetNDF()));
583   parsFit->AddText(Form("Fit :    %4.4f +- %4.4f",fOuterChip->GetParameter(0),fOuterChip->GetParError(0)));
584   ratiochipEvents ->GetListOfFunctions()->Add(parsFit);
585   
586   
587   //------------------------------------------
588   
589   TH1D * diffsClus[2];
590   diffsClus[0]= new TH1D("diffsL1Clus"," ",80,-0.2,0.2);
591   diffsClus[1]= new TH1D("diffsL2Clus"," ",80,-0.2,0.2);
592   for(Int_t ibin=0; ibin<1200; ibin++){
593
594     if(ibin<400){
595       if(ratiochip->GetBinContent(ibin+1)>0) diffsClus[0]->Fill(ratiochip->GetBinContent(ibin+1)-f1innerChip->GetParameter(0));
596     }else {
597
598       if(ratiochip->GetBinContent(ibin+1)>0) diffsClus[1]->Fill(ratiochip->GetBinContent(ibin+1)-f1outerChip->GetParameter(0));
599
600     }
601
602
603   }
604
605
606   TCanvas * pullsClus = new TCanvas("pullsClus","pullsClus",1200,600);
607   pullsClus->Divide(2,1);
608   pullsClus->cd(1);
609   diffsClus[0]->SetTitle("dispersion Layer 1");
610   diffsClus[0]->Rebin(2);
611   diffsClus[0]->Fit("gaus","","",-0.2,0.2);
612   if(diffsClus[0]->GetFunction("gaus")) diffsClus[0]->GetFunction("gaus")->SetLineColor(kBlue);
613   diffsClus[0]->Draw();
614   pullsClus->cd(2);
615   diffsClus[1]->SetTitle("dispersion Layer 2");
616   diffsClus[1]->Rebin(2);
617   diffsClus[1]->Fit("gaus","","",-0.2,0.2);
618   if(diffsClus[1]->GetFunction("gaus")) diffsClus[1]->GetFunction("gaus")->SetLineColor(kBlue);
619   diffsClus[1]->Draw();
620
621 }