]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/pid/doeffKa.C
Change trigger selection logic to include overlaps between J1 and J2 (from Marta)
[u/mrichter/AliRoot.git] / PWGPP / pid / doeffKa.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 LoadLib();
12 void doeffKaUser(Int_t pos,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
13 void doeffKa(Int_t pos=1,Float_t prob=0.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
14 TH2F *GetHistoKap(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);
15 TH2F *GetHistoKan(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);
16 void fit(TH1D *h,Float_t *a=NULL,char *opt="",char *opt2="",Float_t pt=1.5);
17 void AddHisto(TH2F *h1,TH2F *h2,Float_t w);
18 TH2F *GetHistoUser(Int_t pos=1,Float_t pt=1,Float_t ptM=1.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
19 TH2F *GetHistoPiUser(Int_t pos=1,Float_t pt=1,Float_t ptM=1.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
20 TH2F *GetHistoKaUser(Int_t pos=1,Float_t pt=1,Float_t ptM=1.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
21 TH2F *GetHistoPrUser(Int_t pos=1,Float_t pt=1,Float_t ptM=1.1,Float_t etaminkp=-0.8,Float_t etamaxkp=0.8);
22
23 TObject* fContPid1;
24 TObject* fContPid2;
25 TObject* fContUser1;
26 TObject* fContUser2;
27 const Int_t nBinPid = 14; // pt,eta, ptPip, ptPin, PPip, PPin, TOF3sigmaPip, TOF3sigmaPin, isPhiTrue, nsigmaPip, nsigmaPin
28 // 0.985 < mass < 1.045 (60) and 0 < centrality < 100 (10)
29 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*/};
30 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};
31 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};
32
33 TF1 *fsign;
34 TF1 *fall;
35 TF1 *fback;
36
37 Int_t ifunc=0;
38
39 Float_t fitmin = 0.99;
40 Float_t fitmax = 1.045;
41
42 Int_t cmin = 1;// min 1
43 Int_t cmax = 10;// max 10
44
45 Float_t weightS = -0.9;
46
47 Int_t rebinsize = 1;
48
49 Int_t parplotted = 2;
50
51 Bool_t isMC = kFALSE; // don't change this (is set automatically)
52 Bool_t selectTrue = kTRUE; // put it to true to remove background (only for MC)
53 Bool_t keepTrue = kFALSE; // put it to false to fit only background (only for MC)
54
55 Bool_t kGoodMatch = kFALSE; // to check good matching
56
57 Bool_t kSigma2vs3 = kFALSE; // to check good matching
58 Bool_t kSigma2vs3TPC = kFALSE; // to check good matching
59
60 Bool_t require5sigma = kFALSE; // don't touch this flag
61
62 Bool_t bayesVsigma = kFALSE; // only to do checks
63
64 Bool_t kTOFmatch = kFALSE; // for combined PID requires TOF matching
65
66 Bool_t kOverAll = kFALSE;
67 Bool_t kOverAllTOFmatch = kFALSE;
68 Bool_t kOverAll2Sigma = kFALSE;
69
70 TH2F *hmatched;
71 TH2F *htracked;
72
73 Bool_t kLoaded=kFALSE;
74 int LoadLib(){
75   weightS = -0.9;
76
77   require5sigma = kFALSE;
78
79   if(! kLoaded){
80     gSystem->Load("libVMC.so");
81     gSystem->Load("libPhysics.so");
82     gSystem->Load("libTree.so");
83     gSystem->Load("libMinuit.so");
84     gSystem->Load("libSTEERBase.so");
85     gSystem->Load("libANALYSIS.so");
86     gSystem->Load("libAOD.so");
87     gSystem->Load("libESD.so");
88     gSystem->Load("libANALYSIS.so");
89     gSystem->Load("libANALYSISalice.so");
90     gSystem->Load("libCORRFW.so");
91     gSystem->Load("libNetx.so");
92     gSystem->Load("libPWGPPpid.so");
93
94     TFile *f = new TFile("AnalysisResults.root");
95     TList *l = (TList *) f->Get("contPhiBayes1");
96     TList *l2 = (TList *) f->Get("contPhiBayes2");
97
98     if(!(l && l2)) return 0;
99
100     fContPid1 = (AliPIDperfContainer *) l->FindObject("contPID");
101     fContPid2 = (AliPIDperfContainer *) l->FindObject("contPID2");
102     fContUser1 = (AliPIDperfContainer *) l->FindObject("contUserPID");
103     fContUser2 = (AliPIDperfContainer *) l->FindObject("contUserPID2");
104     hmatched = (TH2F *) l2->FindObject("hMatchKa"); 
105     htracked = (TH2F *) l2->FindObject("hTrackingKa"); 
106   }
107   kLoaded = kTRUE;
108
109   // check if MC
110   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]};
111   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]};
112
113   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid1;
114   TH1D *h = tmp->GetQA(0, x, x2)->ProjectionX("checkMC");
115
116   if(h->GetEntries()) isMC = kTRUE;
117   else isMC=kFALSE;
118
119   if(!isMC){
120     selectTrue = kFALSE;
121     keepTrue = kTRUE;
122   }
123   else{
124     printf("MC truth found!!!!!!\nIt is MC!!!!!!");
125   }
126
127   fsign = new TF1("fsign","[0]*TMath::Voigt(x-[1],[3],[2])*(x>0.987)*(x > 1.005 && x < 1.035 || [4])",fitmin,fitmax);
128   fback = new TF1("fback","([0]*sqrt(x-0.987) + [1]*(x-0.987) + [2]*sqrt(x-0.987)*(x-0.987) +[3]*(x-0.987)*(x-0.987)+[4]*(x-0.987)*(x-0.987)*sqrt(x-0.987))*(x>0.987)",fitmin,fitmax);
129   fall = new TF1("fall","([0]*TMath::Voigt(x-[1],[3],[2])*(x > 1.005 && x < 1.035 || [9]) + [4]*sqrt(x-0.987) + [5]*(x-0.987) + [6]*sqrt(x-0.987)*(x-0.987) +[7]*(x-0.987)*(x-0.987)+[8]*(x-0.987)*(x-0.987)*sqrt(x-0.987))*(x>0.987)",0.987,1.05);
130
131   if(isMC){
132     fsign->SetParameter(4,0);
133     fall->FixParameter(9,0);
134   }
135   else{
136     fsign->SetParameter(4,1);
137     fall->FixParameter(9,1);
138   }
139
140   fsign->SetLineColor(2);
141   fback->SetLineColor(4);
142
143   if(kSigma2vs3){
144     kGoodMatch=kFALSE;
145     kOverAll = 0;
146   }
147
148   if(kSigma2vs3TPC){
149     kGoodMatch=kFALSE;
150     kOverAll = 0;
151   }
152
153   if(bayesVsigma){
154     kOverAll = 0;
155     kGoodMatch=kFALSE;
156     kSigma2vs3=kFALSE;
157     kSigma2vs3TPC=kFALSE;
158     kTOFmatch=kTRUE;
159     weightS = -0.5;
160   }
161   if(kOverAll){
162     weightS = -0.5;
163   }
164
165   return 1;
166 }
167 void doeffKaUser(Int_t pos,Float_t etaminkp,Float_t etamaxkp){
168   Int_t nptbin = binPid[2];
169   Float_t minptbin = xmin[2];
170   Float_t maxptbin = xmax[2];
171   
172   TCanvas *c1 = new TCanvas();
173   c1->Divide((nptbin+1)/2,2);
174
175   Double_t xx[50],yyPi[50],yyKa[50],yyPr[50];
176   Double_t exx[50],eyyPi[50],eyyKa[50],eyyPr[50];
177
178   TH2F *hh;
179   TH1D *h;
180
181   Float_t b[100][3];
182   Float_t bPi[100][3];
183   Float_t bKa[100][3];
184   Float_t bPr[100][3];
185
186   char name[100];
187
188   for(Int_t i=0;i < nptbin;i++){
189     c1->cd(i+1);
190     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
191     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
192     
193     xx[i] = (ptmin+ptmax)/2;
194     exx[i] = (-ptmin+ptmax)/2;
195     
196     hh=GetHistoUser(pos,ptmin,ptmax,etaminkp,etamaxkp);
197     sprintf(name,"all%i",i);
198     h = hh->ProjectionX(name,cmin,cmax);
199     Int_t ntrial = 0;
200     Float_t chi2 = 10000;
201     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
202       fit(h,b[i],"WW","",xx[i]);
203       fit(h,b[i],"","",xx[i]);
204       ntrial++;
205       chi2 = b[i][2];
206     }
207     printf("%i) %f +/- %f\n",i,b[i][0],b[i][1]);
208
209     hh=GetHistoPiUser(pos,ptmin,ptmax,etaminkp,etamaxkp);
210     sprintf(name,"pi%i",i);
211     h = hh->ProjectionX(name,cmin,cmax);
212     ntrial = 0;
213     chi2 = 10000;
214     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
215       fit(h,bPi[i],"WW","",xx[i]);
216       fit(h,bPi[i],"","",xx[i]);
217       ntrial++;
218       chi2 = bPi[i][2];
219     }
220     printf("pi) %f +/- %f\n",bPi[i][0],bPi[i][1]);
221
222     hh=GetHistoKaUser(pos,ptmin,ptmax,etaminkp,etamaxkp);
223     sprintf(name,"ka%i",i);
224     h = hh->ProjectionX(name,cmin,cmax);
225     ntrial = 0;
226     chi2 = 10000;
227     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
228       fit(h,bKa[i],"WW","",xx[i]);
229       fit(h,bKa[i],"","",xx[i]);
230       ntrial++;
231       chi2 = bKa[i][2];
232     }
233     printf("ka) %f +/- %f\n",bKa[i][0],bKa[i][1]);
234
235     hh=GetHistoPrUser(pos,ptmin,ptmax,etaminkp,etamaxkp);
236     sprintf(name,"pr%i",i);
237     h = hh->ProjectionX(name,cmin,cmax);
238     ntrial = 0;
239     chi2 = 10000;
240     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
241       fit(h,bPr[i],"WW","",xx[i]);
242       fit(h,bPr[i],"","",xx[i]);
243       ntrial++;
244       chi2 = bPr[i][2];
245     }
246     printf("pr) %f +/- %f\n",bPr[i][0],bPr[i][1]);
247
248     yyPi[i] = bPi[i][0] / b[i][0];
249     yyKa[i] = bKa[i][0] / b[i][0];
250     yyPr[i] = bPr[i][0] / b[i][0];
251
252     eyyPi[i] = bPi[i][1]/bPi[i][0]*yyPi[i];
253     eyyKa[i] = bKa[i][1]/bKa[i][0]*yyKa[i];
254     eyyPr[i] = bPr[i][1]/bPr[i][0]*yyPr[i];
255   }
256
257   /*TCanvas *c2 =*/ new TCanvas();
258   TGraphErrors *gKa = new TGraphErrors(nptbin,xx,yyKa,exx,eyyKa);
259   gKa->Draw("AP");
260   gKa->SetLineColor(1);
261   gKa->SetMarkerColor(1);
262   gKa->SetMarkerStyle(21);
263
264   TGraphErrors *gPi = new TGraphErrors(nptbin,xx,yyPi,exx,eyyPi);
265   gPi->Draw("P");
266   gPi->SetLineColor(4);
267   gPi->SetMarkerColor(4);
268   gPi->SetMarkerStyle(20);
269
270   TGraphErrors *gPr = new TGraphErrors(nptbin,xx,yyPr,exx,eyyPr);
271   gPr->Draw("P");
272   gPr->SetLineColor(2);
273   gPr->SetMarkerColor(2);
274   gPr->SetMarkerStyle(22);
275
276   if(pos) sprintf(name,"phiUserAnalPos_%3.1f-%3.1f_%i-%i.root",etaminkp,etamaxkp,cmin,cmax);
277   else sprintf(name,"phiUserAnalNeg_%3.1f-%3.1f_%i-%i.root",etaminkp,etamaxkp,cmin,cmax);
278
279   gPi->SetName("piSelected");
280   gKa->SetName("kaSelected");
281   gPr->SetName("prSelected");
282
283   TFile *fout = new TFile(name,"RECREATE");
284   gPi->Write();
285   gKa->Write();
286   gPr->Write();
287   fout->Close();
288 }
289
290 TH2F *GetHistoUser(Int_t pos,Float_t pt,Float_t ptM,Float_t etaminkp,Float_t etamaxkp){
291   //  Int_t binUser[nBinUser] = {8/*Eta*/,20/*pt*/,2/*istrue*/,4/*whatSelection*/,1/*DeltaPhi*/,1/*Psi*/};
292   Float_t x[] = {etaminkp+0.0001,pt+0.0001,isMC,0.0001,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
293   Float_t x2[] = {etamaxkp-0.0001,ptM-0.0001,isMC,2.9999,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
294
295   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContUser1;
296   if(!pos) tmp = (AliPIDperfContainer *) fContUser2;
297
298   TH2F *h = tmp->GetQA(0, x, x2);
299
300   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
301   h->GetYaxis()->SetTitle("centrality [%]");
302
303   return h;
304
305 }
306
307
308 TH2F *GetHistoPiUser(Int_t pos,Float_t pt,Float_t ptM,Float_t etaminkp,Float_t etamaxkp){
309   //  Int_t binUser[nBinUser] = {8/*Eta*/,20/*pt*/,2/*istrue*/,4/*whatSelection*/,1/*DeltaPhi*/,1/*Psi*/};
310   Float_t x[] = {etaminkp+0.0001,pt+0.0001,isMC,1.0001,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
311   Float_t x2[] = {etamaxkp-0.0001,ptM-0.0001,isMC,1.0002,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
312
313   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContUser1;
314   if(!pos) tmp = (AliPIDperfContainer *) fContUser2;
315
316   TH2F *h = tmp->GetQA(0, x, x2);
317
318   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
319   h->GetYaxis()->SetTitle("centrality [%]");
320
321   return h;
322
323 }
324
325
326 TH2F *GetHistoKaUser(Int_t pos,Float_t pt,Float_t ptM,Float_t etaminkp,Float_t etamaxkp){
327   //  Int_t binUser[nBinUser] = {8/*Eta*/,20/*pt*/,2/*istrue*/,4/*whatSelection*/,1/*DeltaPhi*/,1/*Psi*/};
328   Float_t x[] = {etaminkp+0.0001,pt+0.0001,isMC,2.0001,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
329   Float_t x2[] = {etamaxkp-0.0001,ptM-0.0001,isMC,2.0002,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
330
331   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContUser1;
332   if(!pos) tmp = (AliPIDperfContainer *) fContUser2;
333
334   TH2F *h = tmp->GetQA(0, x, x2);
335
336   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
337   h->GetYaxis()->SetTitle("centrality [%]");
338
339   return h;
340
341 }
342
343
344 TH2F *GetHistoPrUser(Int_t pos,Float_t pt,Float_t ptM,Float_t etaminkp,Float_t etamaxkp){
345   //  Int_t binUser[nBinUser] = {8/*Eta*/,20/*pt*/,2/*istrue*/,4/*whatSelection*/,1/*DeltaPhi*/,1/*Psi*/};
346   Float_t x[] = {etaminkp+0.0001,pt+0.0001,isMC,3.0001,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
347   Float_t x2[] = {etamaxkp-0.0001,ptM-0.0001,isMC,3.0002,-TMath::Pi(),TMath::Pi(),-TMath::Pi()/2,TMath::Pi()/2};
348
349   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContUser1;
350   if(!pos) tmp = (AliPIDperfContainer *) fContUser2;
351
352   TH2F *h = tmp->GetQA(0, x, x2);
353
354   h->GetXaxis()->SetTitle("M_{K^{0}_{s}} (GeV/#it{c}^{2})");
355   h->GetYaxis()->SetTitle("centrality [%]");
356
357   return h;
358
359 }
360
361
362 void doeffKa(Int_t pos,Float_t prob,Float_t etaminkp,Float_t etamaxkp){
363   LoadLib();
364   TH1D *hm = hmatched->ProjectionX("matchingKaEff",cmin,cmax);
365   TH1D *ht = htracked->ProjectionX("tracking",cmin,cmax);
366
367   hm->GetYaxis()->SetTitle("TOF matching eff.");
368   hm->SetTitle("Using probability as weights");
369
370   hm->Sumw2();
371   ht->Sumw2();
372
373   hm->Divide(hm,ht,1,1,"B");
374
375  
376   Int_t nptbin = binPid[2];
377   Float_t minptbin = xmin[2];
378   Float_t maxptbin = xmax[2];
379
380   if(pos == 0){
381     nptbin = binPid[3];
382     minptbin = xmin[3];
383     maxptbin = xmax[3];
384   }
385
386   if(prob > 0.1999){
387     kGoodMatch = kFALSE;
388     kSigma2vs3 = kFALSE;
389     kSigma2vs3TPC=kFALSE;
390     if(! kOverAll) require5sigma = kTRUE;
391     if(!isMC && !kOverAll) weightS = -0.9;
392   }
393
394   TCanvas *c1 = new TCanvas();
395   c1->Divide((nptbin+1)/2,2);
396   TH2F *hh,*hh2;
397   TH1D *h;
398   char name[100];
399   Float_t b[50][3];
400
401   Double_t xx[50],yy[50];
402   Double_t exx[50],eyy[50];
403
404   for(Int_t i=0;i < nptbin;i++){
405     c1->cd(i+1);//->SetLogy();
406     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
407     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
408
409     xx[i] = (ptmin+ptmax)/2;
410     exx[i] = (-ptmin+ptmax)/2;
411
412     Float_t pp=0.1;
413     if(prob < 0.2) pp = 0.;
414     if(pos) hh=GetHistoKap(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
415     else hh=GetHistoKan(ptmin,ptmax,pp,0.0);
416     sprintf(name,"TOF matched: %f < p_{T} < %f GeV/#it{c}",ptmin,ptmax);
417     hh->SetTitle(name);
418     sprintf(name,"hNoPid%i",i);
419     
420     pp=prob;
421     if(prob < 0.2) pp = 0.1;
422     if(pos) hh2=GetHistoKap(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
423     else hh2=GetHistoKan(ptmin,ptmax,pp,0.0);
424     AddHisto(hh,hh2,weightS);
425
426     h = hh->ProjectionX(name,cmin,cmax);
427     h->RebinX(rebinsize);
428     h->Draw("ERR");
429     h->SetMarkerStyle(24);
430     b[i][0]=-1;
431     Int_t ntrial = 0;
432     Float_t chi2 = 10000;
433     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
434       fit(h,b[i],"WW","",xx[i]);
435       c1->Update();
436 //       getchar();
437       fit(h,b[i],"","",xx[i]);
438       ntrial++;
439       chi2 = b[i][2];
440       printf("chi2 = %f\n",chi2);
441       c1->Update();
442 //       getchar();
443       
444     }
445
446     yy[i] = fall->GetParameter(parplotted);
447     eyy[i] = fall->GetParError(parplotted);
448   }
449
450   TGraphErrors *gpar = new TGraphErrors(nptbin,xx,yy,exx,eyy);
451   c1->cd(8);
452 //   gpar->Draw("AP");
453   gpar->SetMarkerStyle(20);
454
455   TCanvas *c2 = new TCanvas();
456   c2->Divide((nptbin+1)/2,2);
457   Float_t b2[50][3];
458
459   for(Int_t i=0;i < nptbin;i++){
460     c2->cd(i+1);
461     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
462     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
463
464     Float_t pp=prob;
465     if(prob < 0.2) pp = 0.1;
466     if(pos) hh=GetHistoKap(ptmin,ptmax,pp,0.0,etaminkp,etamaxkp);
467     else hh=GetHistoKan(ptmin,ptmax,pp,0.0);
468     sprintf(name,"P_{TOF} > 0.8: %f < p_{T} < %f GeV/#it{c}",ptmin,ptmax);
469     hh->SetTitle(name);
470     sprintf(name,"hPid60_%i",i);
471     h = hh->ProjectionX(name,cmin,cmax);
472     h->RebinX(rebinsize);
473     h->Draw("ERR");
474     h->SetMarkerStyle(24);
475     b2[i][0]=-1;
476     Int_t ntrial = 0;
477     Float_t chi2 = 10000;
478     while(ntrial < 3 && (chi2 > 20 + 1000*selectTrue)){
479       fit(h,b2[i],"WW","");
480       fit(h,b2[i],"","");
481       ntrial++;
482       chi2 = b2[i][2];
483       printf("chi2 = %f\n",chi2);
484     }
485     yy[i] = fall->GetParameter(parplotted);
486     eyy[i] = fall->GetParError(parplotted);
487
488   }
489
490   TGraphErrors *gpar2 = new TGraphErrors(nptbin,xx,yy,exx,eyy);
491   c2->cd(8);
492 //   gpar2->Draw("AP");
493   gpar2->SetMarkerStyle(20);
494   
495   Double_t xpt[50],expt[50],eff[50],efferr[50];
496   for(Int_t i=0;i<nptbin;i++){
497     printf("%f +/- %f -  %f +/- %f\n",b[i][0],b[i][1],b2[i][0],b2[i][1]);
498
499     Float_t ptmin = minptbin+(maxptbin-minptbin)/nptbin*(i);
500     Float_t ptmax = minptbin+(maxptbin-minptbin)/nptbin*(i+1);
501
502     xpt[i] = (ptmin+ptmax)/2;
503     expt[i] = (-ptmin+ptmax)/2;
504     eff[i] = b2[i][0]/(b[i][0]-b2[i][0]*weightS);
505
506     b[i][0] = b[i][0]-b2[i][0]*weightS;
507
508     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]));
509
510     if(TMath::Abs(efferr[i]) > 1)efferr[i]=1;
511   }
512   new TCanvas();
513   TGraphErrors *geff = new TGraphErrors(nptbin,xpt,eff,expt,efferr);
514   geff->Draw("AP");
515
516   char flag[100];
517   flag[0] = '\0';
518
519   if(isMC){
520     if(selectTrue) sprintf(flag,"true");
521     else if(!keepTrue) sprintf(flag,"back");
522   }
523
524   char flag2[100];
525   flag2[0] = '\0';
526
527   Bool_t kWriteME = kFALSE;
528
529   char etarange[100];
530   sprintf(etarange,"_%.1f-%.1f_",etaminkp,etamaxkp);
531
532   if(kGoodMatch)
533     sprintf(flag2,"GM");
534
535   if(bayesVsigma)
536     sprintf(flag2,"BayesVsSigma");
537
538   if(kSigma2vs3)
539     sprintf(flag2,"Sigma2vs3");
540
541   if(kSigma2vs3TPC)
542     sprintf(flag2,"Sigma2vs3TPC");
543
544   if(kOverAll)
545     sprintf(flag2,"OverAll");
546   if(kOverAllTOFmatch)
547     sprintf(flag2,"OverAllTOF"); 
548   if(kOverAll2Sigma)
549     sprintf(flag2,"OverAll2sigma"); 
550
551   if(pos){
552     if(prob >=0.2) sprintf(name,"kaonPos%sP%iEff%i_%i%s%s.root",etarange,Int_t(prob*100),(cmin-1)*10,cmax*10,flag,flag2);
553     else{
554       sprintf(name,"kaonPos%sMatchEff%i_%i%s%s.root",etarange,(cmin-1)*10,cmax*10,flag,flag2);
555       if(!(kOverAll || bayesVsigma || kGoodMatch || kSigma2vs3 || kSigma2vs3TPC)) kWriteME = kTRUE;
556     }
557   }
558   else{
559     if(prob >=0.2) sprintf(name,"kaonNeg%sP%iEff%i_%i%s%s.root",etarange,Int_t(prob*100),(cmin-1)*10,cmax*10,flag,flag2);
560     else sprintf(name,"kaonNeg%sMatchEff%i_%i%s%s.root",etarange,(cmin-1)*10,cmax*10,flag,flag2);
561   }
562
563   geff->SetTitle("K efficiency (from #phi);p_{T} (GeV/#it{c};efficiency");
564   TFile *fout = new TFile(name,"RECREATE");
565   geff->Write();
566   if(kWriteME) hm->Write();
567   fout->Close();
568
569   if(kWriteME) hm->Draw("SAME");
570 }
571
572 TH2F *GetHistoKap(Float_t pt,Float_t ptM,Float_t pMinkp,Float_t pMinkn,Float_t etaminkp,Float_t etamaxkp){
573
574   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]};
575   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]};
576
577   if(kOverAll){
578     x[6] = 0.0001;
579     x2[9] = 4.9;
580   }
581
582   if(kOverAllTOFmatch && pMinkp > 0.19){
583     x[6] = 1.0001;
584   }
585
586   if(kOverAll2Sigma && pMinkp > 0.09){
587     x2[9] = 2;
588     x[6] = 1.0001;
589   }
590
591   if(kGoodMatch){
592     x[6] = 1.0001;
593     if(pMinkp > 0)
594       x2[9] = 4.9;
595       
596   }
597
598   if(kTOFmatch){
599     x[6] = 1.0001;
600   }
601
602   if(kSigma2vs3){
603     x[6] = 1.0001;
604     x2[9] = 3;
605     if(pMinkp > 0)
606       x2[9] = 2;
607   }
608
609   if(kSigma2vs3TPC){
610     x[6] = 0.0001;
611     x2[6] = 0.0002;
612     x2[9] = 3;
613     if(pMinkp > 0)
614       x2[9] = 2;
615   }
616
617   if(bayesVsigma){
618     if(pMinkp > 0){
619       x[4] = 0.2001;
620       x2[9] = 5;
621     }
622     else{
623       x2[9] = 3;
624     }
625
626     
627   }
628
629   if(require5sigma) x2[9] = 4.9;
630
631   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid1;
632
633   TH2F *h = tmp->GetQA(0, x, x2);
634
635   h->GetXaxis()->SetTitle("M_{#phi} (GeV/#it{c}^{2})");
636   h->GetYaxis()->SetTitle("centrality [%]");
637
638   return h;
639 }
640
641 TH2F *GetHistoKan(Float_t pt,Float_t ptM,Float_t pMinkn,Float_t pMinkp,Float_t etaminkp,Float_t etamaxkp){
642
643   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]};
644   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]};
645
646  if(kOverAll){
647     x[7] = 0.0001;
648     x2[10] = 4.9;
649   }
650
651   if(kOverAllTOFmatch && pMinkn > 0.19){
652     x[7] = 1.0001;
653   }
654
655   if(kOverAll2Sigma && pMinkn > 0.09){
656     x2[10] = 2;
657     x[7] = 1.0001;
658   }
659
660   if(kGoodMatch){
661     x[7] = 1.0001;
662     if(pMinkn > 0)
663       x2[10] = 4.9;
664       
665   }
666
667   if(kTOFmatch){
668     x[7] = 1.0001;
669   }
670
671   if(kSigma2vs3){
672     x[7] = 1.0001;
673     x2[10] = 3;
674     if(pMinkn > 0)
675       x2[10] = 2;
676   }
677
678   if(kSigma2vs3TPC){
679     x[7] = 0.0001;
680     x2[7] = 0.0002;
681     x2[10] = 3;
682     if(pMinkn > 0)
683       x2[10] = 2;
684   }
685  
686  if(bayesVsigma){
687     if(pMinkn > 0){
688       x[5] = 0.2001;
689       x2[10] = 5;
690     }
691     else{
692       x2[10] = 3;
693     }    
694   }
695
696   if(require5sigma) x2[10] = 4.9;
697
698   AliPIDperfContainer *tmp = (AliPIDperfContainer *) fContPid2;
699
700   TH2F *h = tmp->GetQA(0, x, x2);
701
702   h->GetXaxis()->SetTitle("M_{#phi} (GeV/#it{c}^{2})");
703   h->GetYaxis()->SetTitle("centrality [%]");
704
705   return h;
706 }
707
708
709 void fit(TH1D *h,Float_t *a,char *opt,char *opt2,Float_t pt){
710   if(h->GetEntries() < 1){
711     if(a){
712       a[0]=0.01;
713       a[1]=1;
714     }
715     return;
716   }
717
718
719  fall->SetParameter(0,100);
720  fall->SetParameter(0,1.01898 + 2.4e-04*pt);
721  fall->SetParameter(2,0.0044);
722  fall->SetParameter(3,0.0015);
723
724  fall->SetParLimits(0,-100,100000);
725  fall->SetParLimits(1,1.01898 + 2.4e-04*pt-1e-03,1.01898 + 2.4e-04*pt+1e-03);
726  fall->SetParLimits(2,0.0005,0.006);
727  fall->SetParLimits(3,0.001,0.0017);
728
729  fall->FixParameter(1,1.01884 + 2.9891e-04*pt);
730  fall->FixParameter(2,0.0044);
731  fall->FixParameter(3,7.57574e-04 + 3.85408e-04*pt);
732
733  fall->ReleaseParameter(4);
734  fall->ReleaseParameter(5);
735  fall->ReleaseParameter(6);
736  fall->ReleaseParameter(7);
737  fall->ReleaseParameter(8);
738
739
740  if(!kGoodMatch && !kSigma2vs3){
741    if(pt > 1.5){
742      fall->FixParameter(7,0);   
743      fall->FixParameter(8,0);   
744    }
745    if(pt > 1.7){
746      fall->FixParameter(6,0);   
747    }
748  }
749
750  if(selectTrue){
751    fall->FixParameter(4,0);
752    fall->FixParameter(5,0);
753    fall->FixParameter(6,0);
754    fall->FixParameter(7,0);
755    fall->FixParameter(8,0);
756  }
757
758  char name[100];
759  TF1 *ftmp=fall;
760
761  TF1 *ftmp2=new TF1(*fsign);
762  sprintf(name,"fsign%i",ifunc);
763  ftmp2->SetName(name);
764
765  TF1 *ftmp3=new TF1(*fback);
766  sprintf(name,"ftmp3%i",ifunc);
767  ftmp3->SetName(name);
768
769  ifunc++;
770
771  h->Fit(ftmp,opt,opt2,fitmin,fitmax);
772  h->Draw("ERR");
773
774  ftmp2->SetParameter(0,ftmp->GetParameter(0));
775  ftmp2->SetParameter(1,ftmp->GetParameter(1));
776  ftmp2->SetParameter(2,ftmp->GetParameter(2));
777  ftmp2->SetParameter(3,ftmp->GetParameter(3));
778  ftmp2->Draw("SAME");
779  ftmp3->SetParameter(0,ftmp->GetParameter(4));
780  ftmp3->SetParameter(1,ftmp->GetParameter(5));
781  ftmp3->SetParameter(2,ftmp->GetParameter(6));
782  ftmp3->SetParameter(3,ftmp->GetParameter(7));
783  ftmp3->SetParameter(4,ftmp->GetParameter(8));
784  ftmp3->Draw("SAME");
785
786  Float_t mean = ftmp->GetParameter(1);
787  Float_t sigma = 0.0044;//TMath::Abs(ftmp->GetParameter(2));
788
789  Float_t signI = ftmp2->Integral(mean-10*sigma,mean+10*sigma)/h->GetBinWidth(1);
790  Float_t backI = ftmp3->Integral(mean-3*sigma,mean+3*sigma)/h->GetBinWidth(1);
791
792  Float_t errI = TMath::Sqrt(ftmp->GetParError(0)*ftmp->GetParError(0)/(0.001+ftmp->GetParameter(0))/(0.001+ftmp->GetParameter(0)));
793
794  printf("signal(5 sigma) = %f +/- %f(fit) +/- %f(stat)\n",signI,errI*signI,TMath::Sqrt(signI));
795  printf("backgr(3sigma) = %f\n",backI);
796  printf("significance(3 sigma) = %f\n",signI/sqrt(signI+backI));
797
798  if(a){
799    a[0]=signI;
800    a[1]=signI*errI*signI*errI + signI;
801    a[1] = TMath::Sqrt(a[1]);
802    if(ftmp->GetNDF()) a[2] = ftmp->GetChisquare()/ftmp->GetNDF();
803
804
805    if(selectTrue){
806      a[0] = h->Integral(1,h->GetNbinsX());
807      a[1] = TMath::Sqrt(a[0]);
808    }
809  }
810 }
811
812 void AddHisto(TH2F *h1,TH2F *h2,Float_t w){
813   Int_t nbinx = h1->GetNbinsX();
814   Int_t nbiny = h1->GetNbinsY();
815
816   for(Int_t i=1;i<=nbinx;i++){
817     for(Int_t j=1;j<=nbiny;j++){
818       Double_t val = h1->GetBinContent(i,j) + h2->GetBinContent(i,j)*w;
819       Float_t err = TMath::Min(TMath::Sqrt(val),val);
820       h1->SetBinContent(i,j,val);
821       h1->SetBinError(i,j,err);
822     }
823   }
824 }