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