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