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