]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/pid/doeffPi.C
improved the macro to extract PID perf (automaticly through a new macro)
[u/mrichter/AliRoot.git] / PWGPP / pid / doeffPi.C
1 #include"TF1.h"
2 #include"TH1D.h"
3 #include"TH2F.h"
4 #include"TMath.h"
5 #include"TSystem.h"
6 #include"TCanvas.h"
7 #include"TFile.h"
8 #include"TGraphErrors.h"
9 #include"AliPIDperfContainer.h"
10
11 Int_t LoadLib();
12 void doeffPi(Int_t pos=1,Float_t prob=0.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
13 TH2F *GetHistoPip(Float_t pt=1,Float_t ptM=1.1,Float_t pMinkp=0,Float_t pMinkn=0.,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
14 TH2F *GetHistoPin(Float_t pt=1,Float_t ptM=1.1,Float_t pMinkn=0,Float_t pMinkp=0.,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
15 void fit(TH1D *h,Float_t *a=NULL,char *opt="",char *opt2="");
16 void AddHisto(TH2F *h1,TH2F *h2,Float_t w);
17
18 TObject* fContPid1;
19 TObject* fContPid2;
20 const Int_t nBinPid = 14; // pt,eta, ptPip, ptPin, PPip, PPin, TOF3sigmaPip, TOF3sigmaPin, isPhiTrue, nsigmaPip, nsigmaPin
21 // 0.985 < mass < 1.045 (60) and 0 < centrality < 100 (10)
22 Int_t binPid[nBinPid] = {1/*ptPhi*/,8/*EtaPi*/,20/*pt+*/,20/*pt-*/,5/*P+*/,1/*P-*/,2/*TOFmatch+*/,2/*TOFmatch-*/,2/*istrue*/,4/*Nsigma+*/,4/*Nsigma-*/,1/*DeltaPhi+*/,1/*DeltaPhi-*/,1/*Psi*/};
23 Float_t xmin[nBinPid] = {1,-0.8,0.3,0.3,0,0,-0.5,-0.5,-0.5,0,0,-TMath::Pi(),-TMath::Pi(),-TMath::Pi()/2};
24 Float_t xmax[nBinPid] = {5,0.8,4.3,4.3,1,1,1.5,1.5,1.5,7.5,7.5,TMath::Pi(),TMath::Pi(),TMath::Pi()/2};
25
26 TF1 *fsign;
27 TF1 *fall;
28 TF1 *fback;
29
30 Int_t ifunc=0;
31
32 Float_t fitmin = 0.3;
33 Float_t fitmax = 0.7;
34
35 Int_t cmin = 1; // min 1
36 Int_t cmax = 10;//max 10
37
38 Float_t weightS = -1.;
39
40 Int_t rebinsize = 1;
41
42 Int_t parplotted = 2;
43
44 Bool_t isMC = kFALSE; // don't change this (is set automatically)
45 Bool_t selectTrue = kTRUE; // put it to true to remove background (only for MC)
46 Bool_t keepTrue = kFALSE; // put it to false to fit only background (only for MC)
47
48 Bool_t kGoodMatch = kFALSE; // to check good matching
49
50 Bool_t kSigma2vs3 = kFALSE; // to check good matching
51
52 Bool_t require5sigma = kFALSE; // don't touch this flag
53
54 Bool_t bayesVsigma = kFALSE; // only to do checks
55
56 Bool_t kTOFmatch = kFALSE; // for combined PID requires TOF matching
57
58 Bool_t kOverAll = kFALSE;
59 Bool_t kOverAllTOFmatch = kFALSE;
60 Bool_t kOverAll2Sigma = kFALSE;
61
62 TH2F *hmatched;
63 TH2F *htracked;
64
65 Bool_t kLoaded=kFALSE;
66
67 Int_t LoadLib(){
68   weightS = -1.;
69
70   if(! kLoaded){
71     gSystem->Load("libVMC.so");
72     gSystem->Load("libPhysics.so");
73     gSystem->Load("libTree.so");
74     gSystem->Load("libMinuit.so");
75     gSystem->Load("libSTEERBase.so");
76     gSystem->Load("libANALYSIS.so");
77     gSystem->Load("libAOD.so");
78     gSystem->Load("libESD.so");
79     gSystem->Load("libANALYSIS.so");
80     gSystem->Load("libANALYSISalice.so");
81     gSystem->Load("libCORRFW.so");
82     gSystem->Load("libNetx.so");
83     gSystem->Load("libPWGPPpid.so");
84
85     TFile *f = new TFile("AnalysisResults.root");
86     TList *l = (TList *) f->Get("contK0sBayes1");
87     TList *l2 = (TList *) f->Get("contK0sBayes2");
88
89     if(! (l && l2)) return 0;
90
91     fContPid1 = (AliPIDperfContainer *) l->FindObject("contPID");
92     fContPid2 = (AliPIDperfContainer *) l->FindObject("contPID2");
93     hmatched = (TH2F *) l2->FindObject("hMatchPi"); 
94     htracked = (TH2F *) l2->FindObject("hTrackingPi"); 
95   }
96   kLoaded = kTRUE;
97
98   // check if MC
99   Float_t x[] = {xmin[0]+0.001,xmin[1]+0.001,xmin[2]+0.001,xmin[3]+0.001,xmin[4]+0.001,xmin[5]+0.001,xmin[6]+0.001,xmin[7]+0.001,1/*trueMC*/,xmin[9],xmin[10]};
100   Float_t x2[] = {xmax[0],xmax[1],xmax[2],xmax[3],xmax[4],xmax[5],xmax[6],xmax[7],xmax[8],xmax[9],xmax[10]};
101
102   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid1;
103   TH1D *h = tmp->GetQA(0, x, x2)->ProjectionX("checkMC");
104
105   if(h->GetEntries()) isMC = kTRUE;
106   else isMC=kFALSE;
107
108   if(!isMC){
109     selectTrue = kFALSE;
110     keepTrue = kTRUE;
111   }
112   else{
113     printf("MC truth found!!!!!!\nIt is MC!!!!!!");
114   }
115
116   fsign = new TF1("fsign","[0]*TMath::Voigt(x-[1],[3],[2])",fitmin,fitmax);
117   fback = new TF1("fback","pol1",fitmin,fitmax);
118   fall = new TF1("fall","[0]*TMath::Voigt(x-[1],[3],[2]) + pol1(4)",fitmin,fitmax);
119
120   fsign->SetLineColor(2);
121   fback->SetLineColor(4);
122
123   if(kSigma2vs3){
124     kGoodMatch=kFALSE;
125     kOverAll = 0;
126   }
127
128   if(bayesVsigma){
129     kOverAll = 0;
130     kGoodMatch=kFALSE;
131     kSigma2vs3=kFALSE;
132     kTOFmatch=kTRUE;
133     weightS = -0.7;
134   }
135
136   return 1;
137 }
138
139 void doeffPi(Int_t pos,Float_t prob,Float_t etaminkp,Float_t etamaxkp){
140   LoadLib();
141   TH1D *hm = hmatched->ProjectionX("matchingPiEff",cmin,cmax);
142   TH1D *ht = htracked->ProjectionX("tracking",cmin,cmax);
143
144   hm->GetYaxis()->SetTitle("TOF matching eff.");
145   hm->SetTitle("Using probability as weights");
146
147   hm->Sumw2();
148   ht->Sumw2();
149
150   hm->Divide(hm,ht,1,1,"B");
151
152   Int_t nptbin = binPid[2];
153   Float_t minptbin = xmin[2];
154   Float_t maxptbin = xmax[2];
155
156   if(pos == 0){
157     nptbin = binPid[3];
158     minptbin = xmin[3];
159     maxptbin = xmax[3];
160   }
161
162   if(prob > 0.1999){
163     kGoodMatch = kFALSE;
164     kSigma2vs3 = kFALSE;
165     if(! kOverAll) require5sigma = kTRUE;
166     if(!isMC) weightS = -0.95;
167   }
168
169   TCanvas *c1 = new TCanvas();
170   c1->Divide((nptbin+1)/2,2);
171   TH2F *hh,*hh2;
172   TH1D *h;
173   char name[100];
174   Float_t b[50][3];
175
176   Double_t xx[50],yy[50];
177   Double_t exx[50],eyy[50];
178
179   for(Int_t i=0;i < nptbin;i++){
180     c1->cd(i+1);//->SetLogy();
181     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
182     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
183
184     xx[i] = (ptmin+ptmax)/2;
185     exx[i] = (-ptmin+ptmax)/2;
186
187     Float_t pp=0.1;
188     if(prob < 0.2) pp = 0.;
189     if(pos) hh=GetHistoPip(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
190     else hh=GetHistoPin(ptmin,ptmax,pp,0.0);
191     sprintf(name,"TOF matched: %f < p_{T} < %f GeV/#it{c}",ptmin,ptmax);
192     hh->SetTitle(name);
193     sprintf(name,"hNoPid%i",i);
194     
195     pp=prob;
196     if(prob < 0.2) pp = 0.1;
197     if(pos) hh2=GetHistoPip(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
198     else hh2=GetHistoPin(ptmin,ptmax,pp,0.0);
199     AddHisto(hh,hh2,weightS);
200
201     h = hh->ProjectionX(name,cmin,cmax);
202     h->RebinX(rebinsize);
203     h->Draw("ERR");
204     h->SetMarkerStyle(24);
205     b[i][0]=-1;
206     Int_t ntrial = 0;
207     Float_t chi2 = 10000;
208     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
209       fit(h,b[i],"WW","");
210       c1->Update();
211 //       getchar();
212       fit(h,b[i],"","");
213       ntrial++;
214       chi2 = b[i][2];
215       printf("chi2 = %f\n",chi2);
216       c1->Update();
217 //       getchar();
218       
219     }
220
221     yy[i] = fall->GetParameter(parplotted);
222     eyy[i] = fall->GetParError(parplotted);
223   }
224
225   TGraphErrors *gpar = new TGraphErrors(nptbin,xx,yy,exx,eyy);
226   c1->cd(8);
227 //   gpar->Draw("AP");
228   gpar->SetMarkerStyle(20);
229
230   TCanvas *c2 = new TCanvas();
231   c2->Divide((nptbin+1)/2,2);
232   Float_t b2[50][3];
233
234   for(Int_t i=0;i < nptbin;i++){
235     c2->cd(i+1);
236     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
237     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
238
239     Float_t pp=prob;
240     if(prob < 0.2) pp = 0.1;
241     if(pos) hh=GetHistoPip(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
242     else hh=GetHistoPin(ptmin,ptmax,pp,0.0);
243     sprintf(name,"P_{TOF} > 0.8: %f < p_{T} < %f GeV/#it{c}",ptmin,ptmax);
244     hh->SetTitle(name);
245     sprintf(name,"hPid60_%i",i);
246     h = hh->ProjectionX(name,cmin,cmax);
247     h->RebinX(rebinsize);
248     h->Draw("ERR");
249     h->SetMarkerStyle(24);
250     b2[i][0]=-1;
251     Int_t ntrial = 0;
252     Float_t chi2 = 10000;
253     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
254       fit(h,b2[i],"WW","");
255       fit(h,b2[i],"","");
256       ntrial++;
257       chi2 = b2[i][2];
258       printf("chi2 = %f\n",chi2);
259     }
260     yy[i] = fall->GetParameter(parplotted);
261     eyy[i] = fall->GetParError(parplotted);
262
263   }
264
265   TGraphErrors *gpar2 = new TGraphErrors(nptbin,xx,yy,exx,eyy);
266   c2->cd(8);
267 //   gpar2->Draw("AP");
268   gpar2->SetMarkerStyle(20);
269   
270   Double_t xpt[50],expt[50],eff[50],efferr[50];
271   for(Int_t i=0;i<nptbin;i++){
272     printf("%f +/- %f -  %f +/- %f\n",b[i][0],b[i][1],b2[i][0],b2[i][1]);
273
274     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
275     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
276
277     xpt[i] = (ptmin+ptmax)/2;
278     expt[i] = (-ptmin+ptmax)/2;
279     eff[i] = b2[i][0]/(b[i][0]-b2[i][0]*weightS);
280
281     b[i][0] = b[i][0]-b2[i][0]*weightS;
282
283     efferr[i] = TMath::Sqrt(b[i][1]*b[i][1]/b[i][0]/b[i][0] + b2[i][1]*b2[i][1]/b2[i][0]/b2[i][0])*(b2[i][0]+b2[i][1])*(1+weightS*(b2[i][0]-b2[i][1])/b[i][0])/b[i][0];//*(1-eff[i]);//der2*der2*(b[i][1]*b[i][1] - b2[i][1]*b2[i][1]));
284
285     if(TMath::Abs(efferr[i]) > 1)efferr[i]=1;
286   }
287   new TCanvas();
288   TGraphErrors *geff = new TGraphErrors(nptbin,xpt,eff,expt,efferr);
289   geff->Draw("AP");
290
291   char flag[100];
292   flag[0] = '\0';
293
294   if(isMC){
295     if(selectTrue) sprintf(flag,"true");
296     else if(!keepTrue) sprintf(flag,"back");
297   }
298
299   Bool_t kWriteME = kFALSE;
300
301   char flag2[100];
302   flag2[0] = '\0';
303
304   char etarange[100];
305   sprintf(etarange,"_%.1f-%.1f_",etaminkp,etamaxkp);
306
307   if(kGoodMatch)
308     sprintf(flag2,"GM");
309
310   if(bayesVsigma)
311     sprintf(flag2,"BayesVsSigma");
312
313   if(kSigma2vs3)
314     sprintf(flag2,"Sigma2vs3");
315
316   if(kOverAll)
317     sprintf(flag2,"OverAll");
318   if(kOverAllTOFmatch)
319     sprintf(flag2,"OverAllTOF");
320   if(kOverAll2Sigma)
321     sprintf(flag2,"OverAll2sigma");
322
323   if(pos){
324     if(prob >=0.2) sprintf(name,"pionPos%sP%iEff%i_%i%s%s.root",etarange,Int_t(prob*100),(cmin-1)*10,cmax*10,flag,flag2);
325     else{
326       sprintf(name,"pionPos%sMatchEff%i_%i%s%s.root",etarange,(cmin-1)*10,cmax*10,flag,flag2);
327       if(!(kOverAll || bayesVsigma || kGoodMatch || kSigma2vs3)) kWriteME = kTRUE;
328     }
329   }
330   else{
331     if(prob >=0.2) sprintf(name,"pionNeg%sP%iEff%i_%i%s%s.root",etarange,Int_t(prob*100),(cmin-1)*10,cmax*10,flag,flag2);
332     else sprintf(name,"pionNeg%sMatchEff%i_%i%s%s.root",etarange,(cmin-1)*10,cmax*10,flag,flag2);
333   }
334
335   geff->SetTitle("#pi efficiency (from K^{0}_{s});p_{T} (GeV/#it{c};efficiency");
336   TFile *fout = new TFile(name,"RECREATE");
337   geff->Write();
338   if(kWriteME) hm->Write();
339   fout->Close();
340
341   if(kWriteME) hm->Draw("SAME");
342 }
343
344 TH2F *GetHistoPip(Float_t pt,Float_t ptM,Float_t pMinkp,Float_t pMinkn,Float_t etaminkp,Float_t etamaxkp){
345
346   Float_t x[] = {xmin[0]+0.001,etaminkp+0.001,pt+0.001,xmin[3]+0.001,pMinkp+0.001,pMinkn+0.001,(pMinkp>0.09)+0.001,kTOFmatch+0.001,selectTrue,xmin[9],xmin[10],xmin[11],xmin[12],xmin[13]};
347   Float_t x2[] = {xmax[0],etamaxkp-0.001,ptM-0.001,xmax[3],xmax[4],xmax[5],xmax[6],xmax[7],keepTrue,xmax[9],xmax[10],xmax[11],xmax[12],xmax[13]};
348
349   if(kOverAll){
350     x[6] = 0.0001;
351     x2[9] = 5.9;
352     if(pMinkp > 0.19) x2[9] = 4.9;
353   }
354
355   if(kOverAllTOFmatch && pMinkp > 0.19){
356     x[6] = 1.0001;
357     x2[9] = 4.9;
358   }
359   
360   if(kOverAll2Sigma && pMinkp > 0.09){
361     x2[9] = 2.;
362     x[6] = 1.0001;
363   }
364
365   if(kGoodMatch){
366     x[6] = 1.0001;
367     if(pMinkp > 0)
368       x2[9] = 4.9;
369       
370   }
371
372   if(kTOFmatch){
373     x[6] = 1.0001;
374   }
375
376   if(kSigma2vs3){
377     x[6] = 1.0001;
378     x2[9] = 3;
379     if(pMinkp > 0)
380       x2[9] = 2;
381   }
382
383   if(bayesVsigma){
384     if(pMinkp > 0){
385       x[4] = 0.2001;
386       x2[9] = 5;
387     }
388     else{
389       x2[9] = 3;
390     }
391
392     
393   }
394
395   if(require5sigma) x2[9] = 4.9;
396
397   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid1;
398
399   TH2F *h = tmp->GetQA(0, x, x2);
400
401   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
402   h->GetYaxis()->SetTitle("centrality [%]");
403
404   return h;
405 }
406
407 TH2F *GetHistoPin(Float_t pt,Float_t ptM,Float_t pMinkn,Float_t pMinkp,Float_t etaminkp,Float_t etamaxkp){
408
409   Float_t x[] = {xmin[0]+0.001,etaminkp+0.001,xmin[2]+0.001,pt+0.001,pMinkp+0.001,pMinkn+0.001,kTOFmatch+0.001,(pMinkn>0.09)+0.001,selectTrue,xmin[9],xmin[10],xmin[11],xmin[12],xmin[13]};
410   Float_t x2[] = {xmax[0],etamaxkp-0.001,xmax[2],ptM-0.001,xmax[4],xmax[5],xmax[6],xmax[7],keepTrue,xmax[9],xmax[10],xmax[11],xmax[12],xmax[13]};
411
412   if(kOverAll){
413     x[7] = 0.0001;
414     x2[10] = 5.9;
415     if(pMinkn > 0.19) x2[10] = 4.9;
416   }
417
418   if(kOverAllTOFmatch && pMinkn > 0.19){
419     x[7] = 1.0001;
420     x2[10] = 4.9;
421   }
422
423   if(kOverAll2Sigma && pMinkn > 0.09){
424     x2[10] = 2;
425     x[7] = 1.0001;
426   }
427
428   if(kGoodMatch){
429     x[7] = 1.0001;
430     if(pMinkn > 0)
431       x2[10] = 4.9;
432       
433   }
434
435   if(kTOFmatch){
436     x[7] = 1.0001;
437   }
438
439   if(kSigma2vs3){
440     x[7] = 1.0001;
441     x2[10] = 3;
442     if(pMinkn > 0)
443       x2[10] = 2;
444   }
445  
446  if(bayesVsigma){
447     if(pMinkn > 0){
448       x[5] = 0.2001;
449       x2[10] = 5;
450     }
451     else{
452       x2[10] = 3;
453     }    
454   }
455
456   if(require5sigma) x2[10] = 4.9;
457
458   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid2;
459
460   TH2F *h = tmp->GetQA(0, x, x2);
461
462   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
463   h->GetYaxis()->SetTitle("centrality [%]");
464
465   return h;
466 }
467
468 void fit(TH1D *h,Float_t *a,char *opt,char *opt2){
469   if(h->Integral(1,h->GetNbinsX()) < 1){
470     if(a){
471       a[0]=0.001;
472       a[1]=1;
473     }
474     return;
475   }
476   
477
478  fall->SetParameter(0,100);
479  fall->SetParameter(1,0.4971);
480  fall->SetParameter(2,0.00006);
481  fall->SetParameter(3,0.0035);
482
483  fall->SetParLimits(0,0.00001,10000);
484  fall->SetParLimits(1,0.4965,0.4985);
485  fall->SetParLimits(2,0.00005,0.001);
486  fall->SetParLimits(3,0.002,0.01);
487
488  fall->ReleaseParameter(4);
489  fall->ReleaseParameter(5);
490
491  if(selectTrue){
492    fall->FixParameter(4,0);
493    fall->FixParameter(5,0);
494  }
495
496  char name[100];
497  TF1 *ftmp=fall;
498
499  TF1 *ftmp2=new TF1(*fsign);
500  sprintf(name,"fsign%i",ifunc);
501  ftmp2->SetName(name);
502
503  TF1 *ftmp3=new TF1(*fback);
504  sprintf(name,"ftmp3%i",ifunc);
505  ftmp3->SetName(name);
506
507  ifunc++;
508
509  h->Fit(ftmp,opt,opt2,fitmin,fitmax);
510  h->Draw("ERR");
511
512  ftmp2->SetParameter(0,ftmp->GetParameter(0));
513  ftmp2->SetParameter(1,ftmp->GetParameter(1));
514  ftmp2->SetParameter(2,ftmp->GetParameter(2));
515  ftmp2->SetParameter(3,ftmp->GetParameter(3));
516  ftmp2->Draw("SAME");
517  ftmp3->SetParameter(0,ftmp->GetParameter(4));
518  ftmp3->SetParameter(1,ftmp->GetParameter(5));
519  ftmp3->Draw("SAME");
520
521  Float_t mean = ftmp->GetParameter(1);
522  Float_t sigma = 0.0044;
523
524  Float_t signI = ftmp2->Integral(mean-10*sigma,mean+10*sigma)/h->GetBinWidth(1);
525  Float_t backI = ftmp3->Integral(mean-3*sigma,mean+3*sigma)/h->GetBinWidth(1);
526
527  Float_t errI = TMath::Sqrt(ftmp->GetParError(0)*ftmp->GetParError(0)/(0.001+ftmp->GetParameter(0))/(0.001+ftmp->GetParameter(0)));
528
529  printf("signal(5 sigma) = %f +/- %f(fit) +/- %f(stat)\n",signI,errI*signI,TMath::Sqrt(signI));
530  printf("backgr(3sigma) = %f\n",backI);
531  printf("significance(3 sigma) = %f\n",signI/sqrt(signI+backI));
532
533  if(a){
534    a[0]=signI;
535    a[1]=signI*errI*signI*errI + signI;
536    a[1] = TMath::Sqrt(a[1]);
537    if(ftmp->GetNDF()) a[2] = ftmp->GetChisquare()/ftmp->GetNDF();
538
539
540    if(selectTrue){
541      a[0] = h->Integral(1,h->GetNbinsX());
542      a[1] = TMath::Sqrt(a[0]);
543    }
544  }
545 }
546
547 void AddHisto(TH2F *h1,TH2F *h2,Float_t w){
548   Int_t nbinx = h1->GetNbinsX();
549   Int_t nbiny = h1->GetNbinsY();
550
551   for(Int_t i=1;i<=nbinx;i++){
552     for(Int_t j=1;j<=nbiny;j++){
553       Double_t val = h1->GetBinContent(i,j) + h2->GetBinContent(i,j)*w;
554       Float_t err = TMath::Min(TMath::Sqrt(val),val);
555       h1->SetBinContent(i,j,val);
556       h1->SetBinError(i,j,err);
557     }
558   }
559 }