Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold7c.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the classes TUnfoldDensity and TUnfoldBinning.
5///
6/// A toy test of the TUnfold package
7///
8///
9/// This example is documented in conference proceedings:
10///
11/// arXiv:1611.01927
12/// 12th Conference on Quark Confinement and the Hadron Spectrum (Confinement XII)
13///
14/// This is an example of unfolding a one-dimensional distribution. It compares
15/// various unfolding methods:
16///
17/// matrix inversion, template fit, L-curve scan,
18/// iterative unfolding, etc
19///
20/// Further details can be found in talk by S.Schmitt at:
21///
22/// XII Quark Confinement and the Hadron Spectrum
23/// 29.8. - 3.9.2016 Thessaloniki, Greece
24/// statictics session (+proceedings)
25///
26/// The example comprises several macros
27/// - testUnfold7a.C create root files with TTree objects for
28/// signal, background and data
29/// - write files testUnfold7_signal.root
30/// testUnfold7_background.root
31/// testUnfold7_data.root
32///
33/// - testUnfold7b.C loop over trees and fill histograms based on the
34/// TUnfoldBinning objects
35/// - read testUnfold7binning.xml
36/// testUnfold7_signal.root
37/// testUnfold7_background.root
38/// testUnfold7_data.root
39///
40/// - write testUnfold7_histograms.root
41///
42/// - testUnfold7c.C run the unfolding
43/// - read testUnfold7_histograms.root
44/// - write testUnfold7_result.root
45/// - write many histograms, to compare various unfolding methods
46///
47/// \macro_output
48/// \macro_image
49/// \macro_code
50///
51/// **Version 17.6, in parallel to changes in TUnfold**
52///
53/// This file is part of TUnfold.
54///
55/// TUnfold is free software: you can redistribute it and/or modify
56/// it under the terms of the GNU General Public License as published by
57/// the Free Software Foundation, either version 3 of the License, or
58/// (at your option) any later version.
59///
60/// TUnfold is distributed in the hope that it will be useful,
61/// but WITHOUT ANY WARRANTY; without even the implied warranty of
62/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
63/// GNU General Public License for more details.
64///
65/// You should have received a copy of the GNU General Public License
66/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
67///
68/// \author Stefan Schmitt DESY, 14.10.2008
69
70#include <iostream>
71#include <cmath>
72#include <map>
73#include <TMath.h>
74#include <TCanvas.h>
75#include <TStyle.h>
76#include <TGraph.h>
77#include <TGraphErrors.h>
78#include <TFile.h>
79#include <TROOT.h>
80#include <TText.h>
81#include <TLine.h>
82#include <TLegend.h>
83#include <TH1.h>
84#include <TF1.h>
85#include <TMatrixD.h>
86#include <TMatrixDSym.h>
87#include <TVectorD.h>
88#include <TMatrixDSymEigen.h>
89#include <TFitResult.h>
90#include <TRandom3.h>
91#include "TUnfoldDensity.h"
92
93using std::vector, std::pair, std::cout;
94
95// #define PRINT_MATRIX_L
96
97#define TEST_INPUT_COVARIANCE
98
99void CreateHistogramCopies(TH1 *h[3],TUnfoldBinning const *binning);
101
102TH2 *AddOverflowXY(TH2 *h,double widthX,double widthY);
103TH1 *AddOverflowX(TH1 *h,double width);
104
105void DrawOverflowX(TH1 *h,double posy);
106void DrawOverflowY(TH1 *h,double posx);
107
108
109double const kLegendFontSize=0.05;
110int kNbinC=0;
111
113void DrawPadEfficiency(TH1 *h);
117 char const *text=nullptr,double tau=0.0,vector<double> const *r=nullptr,
118 TF1 *f=nullptr);
120 vector<pair<TF1*,vector<double> > > const *table);
121
122TFitResultPtr DoFit(TH1 *h,TH2 *rho,TH1 *truth,char const *text,
123 vector<pair<TF1*,vector<double> > > &table,int niter=0);
124
125void GetNiterGraphs(int iFirst,int iLast,vector<pair<TF1*,
126 vector<double> > > const &table,int color,
127 TGraph *graph[4],int style);
128void GetNiterHist(int ifit,vector<pair<TF1*,vector<double> > > const &table,
129 TH1 *hist[4],int color,int style,int fillStyle);
130
131#ifdef WITH_IDS
133
137#endif
138
140
141void testUnfold7c()
142{
144 // switch on histogram errors
146
147 gStyle->SetOptStat(0);
148
149 g_rnd=new TRandom3(4711);
150
151 //==============================================
152 // step 1 : open output file
153 TFile *outputFile=new TFile("testUnfold7_results.root","recreate");
154
155 //==============================================
156 // step 2 : read binning schemes and input histograms
157 TFile *inputFile=new TFile("testUnfold7_histograms.root");
158
159 outputFile->cd();
160
162
163 inputFile->GetObject("fine",fineBinning);
164 inputFile->GetObject("coarse",coarseBinning);
165
166 if((!fineBinning)||(!coarseBinning)) {
167 cout<<"problem to read binning schemes\n";
168 }
169
170 // save binning schemes to output file
171 fineBinning->Write();
172 coarseBinning->Write();
173
174 // read histograms
175#define READ(TYPE,binning,name) \
176 TYPE *name[3]; inputFile->GetObject(#name,name[0]); \
177 name[0]->Write(); \
178 if(!name[0]) cout<<"Error reading " #name "\n"; \
179 CreateHistogramCopies(name,binning);
180
181 outputFile->cd();
182
188
194
197
201 //TH2 *histRhoCLCurve;
202 TH2 *histProbC;
203 double tauMin=1.e-4;
204 double tauMax=1.e-1;
205 double fBgr=1.0; // 0.2/0.25;
206 double biasScale=1.0;
208
209 //double tauC;
210 {
214 mode,
219 tunfoldC->SetInput(histDataRecC[0],biasScale);
220 tunfoldC->SubtractBackground(histMcbgrRecC[0],"BGR",fBgr,0.0);
221 tunfoldC->DoUnfold(0.);
222 histOutputCtau0[0]=tunfoldC->GetOutput("histOutputCtau0");
223 histRhoCtau0=tunfoldC->GetRhoIJtotal("histRhoCtau0");
225 tunfoldC->ScanLcurve(50,tauMin,tauMax,nullptr);
226 /* tauC= */tunfoldC->GetTau();
227 //tunfoldC->ScanTau(50,1.E-7,1.E-1,0,TUnfoldDensity::kEScanTauRhoAvg);
228 histOutputCLCurve[0]=tunfoldC->GetOutput("histOutputCLCurve");
229 /* histRhoCLCurve= */tunfoldC->GetRhoIJtotal("histRhoCLCurve");
231 histProbC=tunfoldC->GetProbabilityMatrix("histProbC",";P_T(gen);P_T(rec)");
232 }
236 //TH2 *histRhoFLCurve;
237 TH2 *histProbF;
238 TGraph *lCurve;
240 tauMin=3.E-4;
241 tauMax=3.E-2;
242 //double tauF;
243 {
247 mode,
252 tunfoldF->SetInput(histDataRecF[0],biasScale);
253 tunfoldF->SubtractBackground(histMcbgrRecF[0],"BGR",fBgr,0.0);
254 tunfoldF->DoUnfold(0.);
255 histOutputFtau0[0]=tunfoldF->GetOutput("histOutputFtau0");
256 histRhoFtau0=tunfoldF->GetRhoIJtotal("histRhoFtau0");
258 tunfoldF->ScanLcurve(50,tauMin,tauMax,nullptr);
259 //tunfoldF->DoUnfold(tauC);
260 /* tauF= */tunfoldF->GetTau();
261 //tunfoldF->ScanTau(50,1.E-7,1.E-1,0,TUnfoldDensity::kEScanTauRhoAvg);
262 histOutputFLCurve[0]=tunfoldF->GetOutput("histOutputFLCurve");
263 /* histRhoFLCurve= */tunfoldF->GetRhoIJtotal("histRhoFLCurve");
265 histProbF=tunfoldF->GetProbabilityMatrix("histProbF",";P_T(gen);P_T(rec)");
266 }
273 TSpline *rhoScan=nullptr;
274 TSpline *logTauCurvature=nullptr;
275
276 double tauFA,tauFArho;
277 {
281 mode,
286 tunfoldFA->SetInput(histDataRecF[0],biasScale);
287 tunfoldFA->SubtractBackground(histMcbgrRecF[0],"BGR",fBgr,0.0);
288 tunfoldFA->DoUnfold(0.);
289 histOutputFAtau0[0]=tunfoldFA->GetOutput("histOutputFAtau0");
290 histRhoFAtau0=tunfoldFA->GetRhoIJtotal("histRhoFAtau0");
293 tauFArho=tunfoldFA->GetTau();
294 histOutputFArho[0]=tunfoldFA->GetOutput("histOutputFArho");
295 histRhoFArho=tunfoldFA->GetRhoIJtotal("histRhoFArho");
297
299 tauFA=tunfoldFA->GetTau();
300 histOutputFALCurve[0]=tunfoldFA->GetOutput("histOutputFALCurve");
301 histRhoFALCurve=tunfoldFA->GetRhoIJtotal("histRhoFALCurve");
303 }
304 lCurve->Write();
305 logTauX->Write();
306 logTauY->Write();
307
308
309 double widthC=coarseBinning->GetBinSize(histProbC->GetNbinsY()+1);
310 double widthF=fineBinning->GetBinSize(histProbF->GetNbinsY()+1);
311
314
315 // efficiency
316 TH1 *histEfficiencyC=histProbCO->ProjectionX("histEfficiencyC");
317 kNbinC=histProbCO->GetNbinsX();
318
319 // reconstructed quantities with overflow (coarse binning)
320 // MC: add signal and bgr
323 histMcbgrRecCO->Scale(fBgr);
324 TH1 *histMcRecCO=(TH1 *)histMcsigRecCO->Clone("histMcRecC0");
327 //TH1 *histDataRecCNO=AddOverflowX(histDataRecC[1],widthC);
328
331 histMcbgrRecFO->Scale(fBgr);
332 TH1 *histMcRecFO=(TH1 *)histMcsigRecFO->Clone("histMcRecF0");
335
336 // truth level with overflow
339
340 // unfolding result with overflow
343 //TH1 *histOutputCLCurveO=AddOverflowX(histOutputCLCurve[2],widthC);
344 //TH2 *histRhoCLCurveO=AddOverflowXY(histRhoCLCurve,widthC,widthC);
349 //TH1 *histOutputFLCurveO=AddOverflowX(histOutputFLCurve[2],widthC);
350 //TH2 *histRhoFLCurveO=AddOverflowXY(histRhoFLCurve,widthC,widthC);
355
356 // bin-by-bin
357 TH2 *histRhoBBBO=(TH2 *)histRhoCtau0O->Clone("histRhoBBBO");
358 for(int i=1;i<=histRhoBBBO->GetNbinsX();i++) {
359 for(int j=1;j<=histRhoBBBO->GetNbinsX();j++) {
360 histRhoBBBO->SetBinContent(i,j,(i==j)?1.:0.);
361 }
362 }
363 TH1 *histDataBgrsub=(TH1 *)histDataRecCO->Clone("histDataBgrsub");
364 histDataBgrsub->Add(histMcbgrRecCO,-fBgr);
365 TH1 *histOutputBBBO=(TH1 *)histDataBgrsub->Clone("histOutputBBBO");
367 histOutputBBBO->Multiply(histMcsigGenO);
368
369 // iterative
370 int niter=1000;
371 cout<<"maximum number of iterations: "<<niter<<"\n";
372
376 histOutputAgo.push_back((TH1*)histMcsigGenO->Clone("histOutputAgo-1"));
377 histOutputAgorep.push_back((TH1*)histMcsigGenO->Clone("histOutputAgorep-1"));
378 histRhoAgo.push_back((TH2*)histRhoBBBO->Clone("histRhoAgo-1"));
379 histRhoAgorep.push_back((TH2*)histRhoBBBO->Clone("histRhoAgorep-1"));
380 nIter.push_back(-1);
381
382 int nx=histProbCO->GetNbinsX();
383 int ny=histProbCO->GetNbinsY();
385 TMatrixD A(ny,nx);
387 for(int i=0;i<nx;i++) {
388 double epsilonI=0.;
389 for(int j=0;j<ny;j++) {
390 epsilonI+= histProbCO->GetBinContent(i+1,j+1);
391 }
392 for(int j=0;j<ny;j++) {
393 double aji=histProbCO->GetBinContent(i+1,j+1);
394 A(j,i)=aji;
396 }
397 }
398 for(int i=0;i<nx;i++) {
400 (histOutputAgo[0]->GetBinError(i+1)
401 *histOutputAgo[0]->GetXaxis()->GetBinWidth(i+1),2.);
402 }
403 for(int i=0;i<ny;i++) {
405 (histDataRecCO->GetBinError(i+1)
406 *histDataRecCO->GetXaxis()->GetBinWidth(i+1),2.);
407 }
408#define NREPLICA 300
413 TVectorD b(nx);
414 for(int nr=0;nr<NREPLICA;nr++) {
415 x[nr]=new TVectorD(nx);
416 y[nr]=new TVectorD(ny);
417 yMb[nr]=new TVectorD(ny);
418 yErr[nr]=new TVectorD(ny);
419 }
420 for(int i=0;i<nx;i++) {
421 (*x[0])(i)=histOutputAgo[0]->GetBinContent(i+1)
422 *histOutputAgo[0]->GetXaxis()->GetBinWidth(i+1);
423 for(int nr=1;nr<NREPLICA;nr++) {
424 (*x[nr])(i)=(*x[0])(i);
425 }
426 }
427 for(int i=0;i<ny;i++) {
428 (*y[0])(i)=histDataRecCO->GetBinContent(i+1)
429 *histDataRecCO->GetXaxis()->GetBinWidth(i+1);
430 for(int nr=1;nr<NREPLICA;nr++) {
431 (*y[nr])(i)=g_rnd->Poisson((*y[0])(i));
432 (*yErr[nr])(i)=TMath::Sqrt((*y[nr])(i));
433 }
434 b(i)=histMcbgrRecCO->GetBinContent(i+1)*
435 histMcbgrRecCO->GetXaxis()->GetBinWidth(i+1);
436 for(int nr=0;nr<NREPLICA;nr++) {
437 (*yMb[nr])(i)=(*y[nr])(i)-b(i);
438 }
439 }
440 for(int iter=0;iter<=niter;iter++) {
441 if(!(iter %100)) cout<<iter<<"\n";
442 for(int nr=0;nr<NREPLICA;nr++) {
443 TVectorD yrec=A*(*x[nr])+b;
445 for(int j=0;j<ny;j++) {
446 yOverYrec(j)=(*y[nr])(j)/yrec(j);
447 }
449 TVectorD xx(nx);
450 for(int i=0;i<nx;i++) {
451 xx(i) = (*x[nr])(i) * f(i);
452 }
453 if(nr==0) {
455 for(int i=0;i<nx;i++) {
456 for(int j=0;j<ny;j++) {
457 xdf_dr(i,j) *= (*x[nr])(i);
458 }
459 }
461 for(int j=0;j<ny;j++) {
462 dr_dxdy(j,nx+j)=1.0/yrec(j);
463 for(int i=0;i<nx;i++) {
464 dr_dxdy(j,i)= -yOverYrec(j)/yrec(j)*A(j,i);
465 }
466 }
468 dxy_dxy.SetSub(0,0,xdf_dr*dr_dxdy);
469 for(int i=0;i<nx;i++) {
470 dxy_dxy(i,i) +=f(i);
471 }
472 for(int i=0;i<ny;i++) {
473 dxy_dxy(nx+i,nx+i) +=1.0;
474 }
477 }
478 (*x[nr])=xx;
479 }
480 if((iter<=25)||
481 ((iter<=100)&&(iter %5==0))||
482 ((iter<=1000)&&(iter %50==0))||
483 (iter %1000==0)) {
484 nIter.push_back(iter);
485 TH1 * h=(TH1*)histOutputAgo[0]->Clone
486 (TString::Format("histOutputAgo%d",iter));
487 histOutputAgo.push_back(h);
488 for(int i=0;i<nx;i++) {
489 double bw=h->GetXaxis()->GetBinWidth(i+1);
490 h->SetBinContent(i+1,(*x[0])(i)/bw);
491 h->SetBinError(i+1,TMath::Sqrt(covAgo(i,i))/bw);
492 }
494 (TString::Format("histRhoAgo%d",iter));
495 histRhoAgo.push_back(h2);
496 for(int i=0;i<nx;i++) {
497 for(int j=0;j<nx;j++) {
498 double rho= covAgo(i,j)/TMath::Sqrt(covAgo(i,i)*covAgo(j,j));
499 if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
500 cout<<"bad error matrix: iter="<<iter<<"\n";
501 exit(0);
502 }
503 h2->SetBinContent(i+1,j+1,rho);
504 }
505 }
506 // error and correlations from replica analysis
507 h=(TH1*)histOutputAgo[0]->Clone
508 (TString::Format("histOutputAgorep%d",iter));
510 (TString::Format("histRhoAgorep%d",iter));
511 histOutputAgorep.push_back(h);
512 histRhoAgorep.push_back(h2);
513
514 TVectorD mean(nx);
515 double w=1./(NREPLICA-1.);
516 for(int nr=1;nr<NREPLICA;nr++) {
517 mean += w* *x[nr];
518 }
520 for(int nr=1;nr<NREPLICA;nr++) {
521 //TMatrixD dx= (*x)-mean;
522 TMatrixD dx(nx,1);
523 for(int i=0;i<nx;i++) {
524 dx(i,0)= (*x[nr])(i)-(*x[0])(i);
525 }
527 }
528
529 for(int i=0;i<nx;i++) {
530 double bw=h->GetXaxis()->GetBinWidth(i+1);
531 h->SetBinContent(i+1,(*x[0])(i)/bw);
532 h->SetBinError(i+1,TMath::Sqrt(covAgorep(i,i))/bw);
533 // cout<<i<<" "<<(*x[0])(i)/bw<<" +/-"<<TMath::Sqrt(covAgorep(i,i))/bw<<" "<<TMath::Sqrt(covAgo(i,i))/bw<<"\n";
534 }
535 for(int i=0;i<nx;i++) {
536 for(int j=0;j<nx;j++) {
537 double rho= covAgorep(i,j)/
539 if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
540 cout<<"bad error matrix: iter="<<iter<<"\n";
541 exit(0);
542 }
543 h2->SetBinContent(i+1,j+1,rho);
544 }
545 }
546 }
547 }
548
549#ifdef WITH_IDS
550 // IDS Malaescu
551 int niterIDS=100;
553 cout<<"IDS number of iterations: "<<niterIDS<<"\n";
556 for(int nr=0;nr<NREPLICA;nr++) {
557 Am_IDS[nr]=new TMatrixD(ny,nx);
558 }
559 for(int iy=0;iy<ny;iy++) {
560 for(int ix=0;ix<nx;ix++) {
561 A_IDS(iy,ix)=histMcsigGenRecC[0]->GetBinContent(ix+1,iy+1);
562 }
563 }
564 double lambdaL=0.;
565 Double_t lambdaUmin = 1.0000002;
566 Double_t lambdaMmin = 0.0000001;
567 Double_t lambdaS = 0.000001;
568 double lambdaU=lambdaUmin;
569 double lambdaM=lambdaMmin;
572 histOutputIDS.push_back((TH1*)histOutputAgo[0]->Clone("histOutputIDS-1"));
573 histRhoIDS.push_back((TH2*)histRhoAgo[0]->Clone("histRhoIDS-1"));
574 histOutputIDS.push_back((TH1*)histOutputAgo[0]->Clone("histOutputIDS0"));
575 histRhoIDS.push_back((TH2*)histRhoAgo[0]->Clone("histRhoIDS0"));
576 for(int iter=1;iter<=niterIDS;iter++) {
577 if(!(iter %10)) cout<<iter<<"\n";
578
579 for(int nr=0;nr<NREPLICA;nr++) {
580 if(iter==1) {
582 } else {
586 }
587 }
588 unsigned ix;
589 for(ix=0;ix<nIter.size();ix++) {
590 if(nIter[ix]==iter) break;
591 }
592 if(ix<nIter.size()) {
593 TH1 * h=(TH1*)histOutputIDS[0]->Clone
594 (TString::Format("histOutputIDS%d",iter));
596 (TString::Format("histRhoIDS%d",iter));
597 histOutputIDS.push_back(h);
598 histRhoIDS.push_back(h2);
599 TVectorD mean(nx);
600 double w=1./(NREPLICA-1.);
601 for(int nr=1;nr<NREPLICA;nr++) {
602 mean += w* (*unfresIDS[nr]);
603 }
605 for(int nr=1;nr<NREPLICA;nr++) {
606 //TMatrixD dx= (*x)-mean;
607 TMatrixD dx(nx,1);
608 for(int i=0;i<nx;i++) {
609 dx(i,0)= (*unfresIDS[nr])(i)-(*unfresIDS[0])(i);
610 }
612 }
613 for(int i=0;i<nx;i++) {
614 double bw=h->GetXaxis()->GetBinWidth(i+1);
615 h->SetBinContent(i+1,(*unfresIDS[0])(i)/bw/
616 histEfficiencyC->GetBinContent(i+1));
617 h->SetBinError(i+1,TMath::Sqrt(covIDSrep(i,i))/bw/
618 histEfficiencyC->GetBinContent(i+1));
619 // cout<<i<<" "<<(*x[0])(i)/bw<<" +/-"<<TMath::Sqrt(covAgorep(i,i))/bw<<" "<<TMath::Sqrt(covAgo(i,i))/bw<<"\n";
620 }
621 for(int i=0;i<nx;i++) {
622 for(int j=0;j<nx;j++) {
623 double rho= covIDSrep(i,j)/
625 if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
626 cout<<"bad error matrix: iter="<<iter<<"\n";
627 exit(0);
628 }
629 h2->SetBinContent(i+1,j+1,rho);
630 }
631 }
632 }
633 }
634#endif
635
636 //double NEdSmc=histDataBgrsub->GetSumOfWeights();
637
639
640 TCanvas *c1=new TCanvas("c1","",600,600);
641 TCanvas *c2sq=new TCanvas("c2sq","",600,600);
642 c2sq->Divide(1,2);
643 TCanvas *c2w=new TCanvas("c2w","",600,300);
644 c2w->Divide(2,1);
645 TCanvas *c4=new TCanvas("c4","",600,600);
646 c4->Divide(2,2);
647 //TCanvas *c3n=new TCanvas("c3n","",600,600);
648 TPad *subn[3];
649 //gROOT->SetStyle("xTimes2");
650 subn[0]= new TPad("subn0","",0.,0.5,1.,1.);
651 //gROOT->SetStyle("square");
652 subn[1]= new TPad("subn1","",0.,0.,0.5,0.5);
653 subn[2]= new TPad("subn2","",0.5,0.0,1.,0.5);
654 for(int i=0;i<3;i++) {
655 subn[i]->SetFillStyle(0);
656 subn[i]->Draw();
657 }
658 TCanvas *c3c=new TCanvas("c3c","",600,600);
659 TPad *subc[3];
660 //gROOT->SetStyle("xTimes2");
661 subc[0]= new TPad("sub0","",0.,0.5,1.,1.);
662 //gROOT->SetStyle("squareCOLZ");
663 subc[1]= new TPad("sub1","",0.,0.,0.5,0.5);
664 //gROOT->SetStyle("square");
665 subc[2]= new TPad("sub2","",0.5,0.0,1.,0.5);
666 for(int i=0;i<3;i++) {
667 subc[i]->SetFillStyle(0);
668 subc[i]->Draw();
669 }
670
671 //=========================== example ==================================
672
673 c2w->cd(1);
675 c2w->cd(2);
676 DrawPadReco(histMcRecCO,histMcbgrRecCO,histDataRecCO,nullptr,nullptr,nullptr);
677 c2w->SaveAs("exampleTR.eps");
678
679 //=========================== example ==================================
680
681 c2w->cd(1);
683 c2w->cd(2);
685 c2w->SaveAs("exampleAE.eps");
686
687 int iFitInversion=table.size();
688 DoFit(histOutputCtau0O,histRhoCtau0O,histDataGenO,"inversion",table);
689 //=========================== inversion ==================================
690
691 subc[0]->cd();
693 &table[table.size()-1].second);
694 subc[1]->cd();
696 subc[2]->cd();
699 c3c->SaveAs("inversion.eps");
700
701
702 DoFit(histOutputFtau0O,histRhoFtau0O,histDataGenO,"template",table);
703 //=========================== template ==================================
704
705 subc[0]->cd();
707 &table[table.size()-1].second);
708 subc[1]->cd();
710 subc[2]->cd();
713 c3c->SaveAs("template.eps");
714
715 DoFit(histOutputFAtau0O,histRhoFAtau0O,histDataGenO,"template+area",table);
716 //=========================== template+area ==================================
717
718 subc[0]->cd();
720 &table[table.size()-1].second);
721 subc[1]->cd();
723 subc[2]->cd();
726 c3c->SaveAs("templateA.eps");
727
728 int iFitFALCurve=table.size();
729 DoFit(histOutputFALCurveO,histRhoFALCurveO,histDataGenO,"Tikhonov+area",table);
730 //=========================== template+area+tikhonov =====================
731 subc[0]->cd();
733 &table[table.size()-1].second);
734 subc[1]->cd();
736 subc[2]->cd();
739 c3c->SaveAs("lcurveFA.eps");
740
741 int iFitFArho=table.size();
742 DoFit(histOutputFArhoO,histRhoFArhoO,histDataGenO,"min(rhomax)",table);
743 //=========================== template+area+tikhonov =====================
744 subc[0]->cd();
746 &table[table.size()-1].second);
747 subc[1]->cd();
749 subc[2]->cd();
752 c3c->SaveAs("rhoscanFA.eps");
753
754 int iFitBinByBin=table.size();
755 DoFit(histOutputBBBO,histRhoBBBO,histDataGenO,"bin-by-bin",table);
756 //=========================== bin-by-bin =================================
757 //c->cd(1);
758 //DrawPadProbability(histProbCO);
759 //c->cd(2);
760 //DrawPadCorrelations(histRhoBBBO,&table);
761 c2sq->cd(1);
764 c2sq->cd(2);
766 &table[table.size()-1].second);
767 c2sq->SaveAs("binbybin.eps");
768
769
770 //=========================== iterative ===================================
771 int iAgoFirstFit=table.size();
772 for(size_t i=1;i<histRhoAgorep.size();i++) {
773 int n=nIter[i];
774 bool isFitted=false;
776 TString::Format("iterative, N=%d",n),table,n);
777 isFitted=true;
778 subc[0]->cd();
780 TString::Format("iterative N=%d",nIter[i]),0.,
781 isFitted ? &table[table.size()-1].second : nullptr);
782 subc[1]->cd();
784 subc[2]->cd();
787 c3c->SaveAs(TString::Format("iterative%d.eps",nIter[i]));
788 }
789 int iAgoLastFit=table.size();
790
791#ifdef WITH_IDS
792 int iIDSFirstFit=table.size();
793
794 //=========================== IDS ===================================
795
796 for(size_t i=2;i<histRhoIDS.size();i++) {
797 int n=nIter[i];
798 bool isFitted=false;
800 TString::Format("IDS, N=%d",n),table,n);
801 isFitted=true;
802 subc[0]->cd();
804 TString::Format("IDS N=%d",nIter[i]),0.,
805 isFitted ? &table[table.size()-1].second : 0);
806 subc[1]->cd();
808 subc[2]->cd();
811 c3c->SaveAs(TString::Format("ids%d.eps",nIter[i]));
812 }
813 int iIDSLastFit=table.size();
814#endif
815
816 int nfit=table.size();
817 TH1D *fitChindf=new TH1D("fitChindf",";algorithm;#chi^{2}/NDF",nfit,0,nfit);
818 TH1D *fitNorm=new TH1D("fitNorm",";algorithm;Landau amplitude [1/GeV]",nfit,0,nfit);
819 TH1D *fitMu=new TH1D("fitMu",";algorithm;Landau #mu [GeV]",nfit,0,nfit);
820 TH1D *fitSigma=new TH1D("fitSigma",";algorithm;Landau #sigma [GeV]",nfit,0,nfit);
821 for(int fit=0;fit<nfit;fit++) {
822 TF1 *f=table[fit].first;
823 vector<double> const &r=table[fit].second;
824 fitChindf->GetXaxis()->SetBinLabel(fit+1,f->GetName());
825 fitNorm->GetXaxis()->SetBinLabel(fit+1,f->GetName());
826 fitMu->GetXaxis()->SetBinLabel(fit+1,f->GetName());
827 fitSigma->GetXaxis()->SetBinLabel(fit+1,f->GetName());
828 double chi2=r[0];
829 double ndf=r[1];
830 fitChindf->SetBinContent(fit+1,chi2/ndf);
831 fitChindf->SetBinError(fit+1,TMath::Sqrt(2./ndf));
832 fitNorm->SetBinContent(fit+1,f->GetParameter(0));
833 fitNorm->SetBinError(fit+1,f->GetParError(0));
834 fitMu->SetBinContent(fit+1,f->GetParameter(1));
835 fitMu->SetBinError(fit+1,f->GetParError(1));
836 fitSigma->SetBinContent(fit+1,f->GetParameter(2));
837 fitSigma->SetBinError(fit+1,f->GetParError(2));
838 cout<<"\""<<f->GetName()<<"\","<<r[2]/r[3]<<","<<r[3]
839 <<","<<TMath::Prob(r[2],r[3])
840 <<","<<chi2/ndf
841 <<","<<ndf
842 <<","<<TMath::Prob(r[0],r[1])
843 <<","<<r[5];
844 for(int i=1;i<3;i++) {
845 cout<<","<<f->GetParameter(i)<<",\"\302\261\","<<f->GetParError(i);
846 }
847 cout<<"\n";
848 }
849
850 //=========================== L-curve ==========================
851 c4->cd(1);
852 lCurve->SetTitle("L curve;log_{10} L_{x};log_{10} L_{y}");
853 lCurve->SetLineColor(kRed);
854 lCurve->Draw("AL");
855 c4->cd(2);
856 gPad->Clear();
857 c4->cd(3);
858 logTauX->SetTitle(";log_{10} #tau;log_{10} L_{x}");
859 logTauX->SetLineColor(kBlue);
860 logTauX->Draw();
861 c4->cd(4);
862 logTauY->SetTitle(";log_{10} #tau;log_{10} L_{y}");
863 logTauY->SetLineColor(kBlue);
864 logTauY->Draw();
865 c4->SaveAs("lcurveL.eps");
866
867 //========================= rho and L-curve scan ===============
868 c4->cd(1);
869 logTauCurvature->SetTitle(";log_{10}(#tau);L curve curvature");
870 logTauCurvature->SetLineColor(kRed);
871 logTauCurvature->Draw();
872 c4->cd(2);
873 rhoScan->SetTitle(";log_{10}(#tau);average(#rho_{i})");
874 rhoScan->SetLineColor(kRed);
875 rhoScan->Draw();
876 c4->cd(3);
878 &table[iFitFALCurve].second);
879 c4->cd(4);
881 &table[iFitFArho].second);
882
883 c4->SaveAs("scanTau.eps");
884
885
890#ifdef WITH_IDS
893#endif
902
903 histNiterInversion[0]->GetYaxis()->SetRangeUser(0.3,500.);
904 histNiterInversion[1]->GetYaxis()->SetRangeUser(-0.1,0.9);
905 histNiterInversion[2]->GetYaxis()->SetRangeUser(5.6,6.3);
906 histNiterInversion[3]->GetYaxis()->SetRangeUser(1.6,2.4);
907
908 TLine *line=nullptr;
909 c1->cd();
910 for(int i=0;i<2;i++) {
911 gPad->Clear();
912 gPad->SetLogx();
913 gPad->SetLogy((i<1));
914 if(! histNiterInversion[i]) continue;
915 histNiterInversion[i]->Draw("][");
916 histNiterFALCurve[i]->Draw("SAME ][");
917 histNiterFArho[i]->Draw("SAME ][");
918 histNiterBinByBin[i]->Draw("SAME ][");
919 graphNiterAgo[i]->Draw("LP");
920 graphNiterAgoBay[i]->SetMarkerStyle(20);
921 graphNiterAgoBay[i]->Draw("P");
922#ifdef WITH_IDS
923 graphNiterIDS[i]->Draw("LP");
924#endif
926 if(i==1) {
927 legend=new TLegend(0.48,0.28,0.87,0.63);
928 } else {
929 legend=new TLegend(0.45,0.5,0.88,0.88);
930 }
931 legend->SetBorderSize(0);
932 legend->SetFillStyle(1001);
933 legend->SetFillColor(kWhite);
934 legend->SetTextSize(kLegendFontSize*0.7);
935 legend->AddEntry( histNiterInversion[0],"inversion","l");
936 legend->AddEntry( histNiterFALCurve[0],"Tikhonov L-curve","l");
937 legend->AddEntry( histNiterFArho[0],"Tikhonov global cor.","l");
938 legend->AddEntry( histNiterBinByBin[0],"bin-by-bin","l");
939 legend->AddEntry( graphNiterAgoBay[0],"\"Bayesian\"","p");
940 legend->AddEntry( graphNiterAgo[0],"iterative","p");
941#ifdef WITH_IDS
942 legend->AddEntry( graphNiterIDS[0],"IDS","p");
943#endif
944 legend->Draw();
945
946 c1->SaveAs(TString::Format("niter%d.eps",i));
947 }
948
949 //c4->cd(1);
950 //DrawPadCorrelations(histRhoFALCurveO,&table);
951 c2sq->cd(1);
953 &table[iFitFALCurve].second,table[iFitFALCurve].first);
954
955 c2sq->cd(2);
956 gPad->SetLogx();
957 gPad->SetLogy(false);
958 histNiterInversion[3]->DrawClone("E2");
959 histNiterInversion[3]->SetFillStyle(0);
960 histNiterInversion[3]->Draw("SAME HIST ][");
961 histNiterFALCurve[3]->DrawClone("SAME E2");
962 histNiterFALCurve[3]->SetFillStyle(0);
963 histNiterFALCurve[3]->Draw("SAME HIST ][");
964 histNiterFArho[3]->DrawClone("SAME E2");
965 histNiterFArho[3]->SetFillStyle(0);
966 histNiterFArho[3]->Draw("SAME HIST ][");
967 histNiterBinByBin[3]->DrawClone("SAME E2");
968 histNiterBinByBin[3]->SetFillStyle(0);
969 histNiterBinByBin[3]->Draw("SAME HIST ][");
970 double yTrue=1.8;
971 line=new TLine(histNiterInversion[3]->GetXaxis()->GetXmin(),
972 yTrue,
973 histNiterInversion[3]->GetXaxis()->GetXmax(),
974 yTrue);
975 line->SetLineWidth(3);
976 line->Draw();
977
978 graphNiterAgo[3]->Draw("LP");
979 graphNiterAgoBay[3]->SetMarkerStyle(20);
980 graphNiterAgoBay[3]->Draw("P");
981#ifdef WITH_IDS
982 graphNiterIDS[3]->Draw("LP");
983#endif
984
986 legend=new TLegend(0.55,0.53,0.95,0.97);
987 legend->SetBorderSize(0);
988 legend->SetFillStyle(1001);
989 legend->SetFillColor(kWhite);
990 legend->SetTextSize(kLegendFontSize);
991 legend->AddEntry( line,"truth","l");
992 legend->AddEntry( histNiterInversion[3],"inversion","l");
993 legend->AddEntry( histNiterFALCurve[3],"Tikhonov L-curve","l");
994 legend->AddEntry( histNiterFArho[3],"Tikhonov global cor.","l");
995 legend->AddEntry( histNiterBinByBin[3],"bin-by-bin","l");
996 legend->AddEntry( graphNiterAgoBay[3],"\"Bayesian\"","p");
997 legend->AddEntry( graphNiterAgo[3],"iterative","p");
998#ifdef WITH_IDS
999 legend->AddEntry( graphNiterIDS[3],"IDS","p");
1000#endif
1001 legend->Draw();
1002
1003 c2sq->SaveAs("fitSigma.eps");
1004
1005 outputFile->Write();
1006
1007 delete outputFile;
1008
1009}
1010
1011void GetNiterGraphs(int iFirst,int iLast,vector<pair<TF1*,
1012 vector<double> > > const &table,int color,
1013 TGraph *graph[4],int style) {
1018 TVectorD mean(iLast-iFirst);
1022 for(int ifit=iFirst;ifit<iLast;ifit++) {
1023 vector<double> const &r=table[ifit].second;
1024 niter(ifit-iFirst)=r[4];
1025 chi2(ifit-iFirst)=r[2]/r[3];
1026 gcor(ifit-iFirst)=r[5];
1027 TF1 const *f=table[ifit].first;
1028 mean(ifit-iFirst)=f->GetParameter(1);
1029 emean(ifit-iFirst)=f->GetParError(1);
1030 sigma(ifit-iFirst)=f->GetParameter(2);
1031 esigma(ifit-iFirst)=f->GetParError(2);
1032 }
1033 graph[0]=new TGraph(niter,chi2);
1034 graph[1]=new TGraph(niter,gcor);
1035 graph[2]=new TGraphErrors(niter,mean,eniter,emean);
1036 graph[3]=new TGraphErrors(niter,sigma,eniter,esigma);
1037 for(int g=0;g<4;g++) {
1038 if(graph[g]) {
1039 graph[g]->SetLineColor(color);
1040 graph[g]->SetMarkerColor(color);
1041 graph[g]->SetMarkerStyle(style);
1042 }
1043 }
1044}
1045
1046void GetNiterHist(int ifit,vector<pair<TF1*,vector<double> > > const &table,
1047 TH1 *hist[4],int color,int style,int fillStyle) {
1048 vector<double> const &r=table[ifit].second;
1049 TF1 const *f=table[ifit].first;
1050 hist[0]=new TH1D(table[ifit].first->GetName()+TString("_chi2"),
1051 ";iteration;unfold-truth #chi^{2}/N_{D.F.}",1,0.2,1500.);
1052 hist[0]->SetBinContent(1,r[2]/r[3]);
1053
1054 hist[1]=new TH1D(table[ifit].first->GetName()+TString("_gcor"),
1055 ";iteration;avg(#rho_{i})",1,0.2,1500.);
1056 hist[1]->SetBinContent(1,r[5]);
1057 hist[2]=new TH1D(table[ifit].first->GetName()+TString("_mu"),
1058 ";iteration;parameter #mu",1,0.2,1500.);
1059 hist[2]->SetBinContent(1,f->GetParameter(1));
1060 hist[2]->SetBinError(1,f->GetParError(1));
1061 hist[3]=new TH1D(table[ifit].first->GetName()+TString("_sigma"),
1062 ";iteration;parameter #sigma",1,0.2,1500.);
1063 hist[3]->SetBinContent(1,f->GetParameter(2));
1064 hist[3]->SetBinError(1,f->GetParError(2));
1065 for(int h=0;h<4;h++) {
1066 if(hist[h]) {
1067 hist[h]->SetLineColor(color);
1068 hist[h]->SetLineStyle(style);
1069 if( hist[h]->GetBinError(1)>0.0) {
1070 hist[h]->SetFillColor(color-10);
1071 hist[h]->SetFillStyle(fillStyle);
1072 }
1073 hist[h]->SetMarkerStyle(0);
1074 }
1075 }
1076}
1077
1078void CreateHistogramCopies(TH1 *h[3],TUnfoldBinning const *binning) {
1079 TString baseName(h[0]->GetName());
1080 Int_t *binMap;
1081 h[1]=binning->CreateHistogram(baseName+"_axis",kTRUE,&binMap);
1082 h[2]=(TH1 *)h[1]->Clone(baseName+"_binw");
1083 Int_t nMax=binning->GetEndBin()+1;
1084 for(Int_t iSrc=0;iSrc<nMax;iSrc++) {
1086 double c=h[0]->GetBinContent(iSrc)+h[1]->GetBinContent(iDest);
1087 double e=TMath::Hypot(h[0]->GetBinError(iSrc),h[1]->GetBinError(iDest));
1088 h[1]->SetBinContent(iDest,c);
1089 h[1]->SetBinError(iDest,e);
1090 h[2]->SetBinContent(iDest,c);
1091 h[2]->SetBinError(iDest,e);
1092 }
1093 for(int iDest=0;iDest<=h[2]->GetNbinsX()+1;iDest++) {
1094 double c=h[2]->GetBinContent(iDest);
1095 double e=h[2]->GetBinError(iDest);
1096 double bw=binning->GetBinSize(iDest);
1097 /* if(bw!=h[2]->GetBinWidth(iDest)) {
1098 cout<<"bin "<<iDest<<" width="<<bw<<" "<<h[2]->GetBinWidth(iDest)
1099 <<"\n";
1100 } */
1101 if(bw>0.0) {
1102 h[2]->SetBinContent(iDest,c/bw);
1103 h[2]->SetBinError(iDest,e/bw);
1104 } else {
1105 }
1106 }
1107}
1108
1110 h[1]=nullptr;
1111 h[2]=nullptr;
1112}
1113
1114TH2 *AddOverflowXY(TH2 *h,double widthX,double widthY) {
1115 // add overflow bin to X-axis
1116 int nx=h->GetNbinsX();
1117 int ny=h->GetNbinsY();
1118 double *xBins=new double[nx+2];
1119 double *yBins=new double[ny+2];
1120 for(int i=1;i<=nx;i++) {
1121 xBins[i-1]=h->GetXaxis()->GetBinLowEdge(i);
1122 }
1123 xBins[nx]=h->GetXaxis()->GetBinUpEdge(nx);
1124 xBins[nx+1]=xBins[nx]+widthX;
1125 for(int i=1;i<=ny;i++) {
1126 yBins[i-1]=h->GetYaxis()->GetBinLowEdge(i);
1127 }
1128 yBins[ny]=h->GetYaxis()->GetBinUpEdge(ny);
1129 yBins[ny+1]=yBins[ny]+widthY;
1130 TString name(h->GetName());
1131 name+="U";
1132 TH2 *r=new TH2D(name,h->GetTitle(),nx+1,xBins,ny+1,yBins);
1133 for(int ix=0;ix<=nx+1;ix++) {
1134 for(int iy=0;iy<=ny+1;iy++) {
1135 r->SetBinContent(ix,iy,h->GetBinContent(ix,iy));
1136 r->SetBinError(ix,iy,h->GetBinError(ix,iy));
1137 }
1138 }
1139 delete [] yBins;
1140 delete [] xBins;
1141 return r;
1142}
1143
1144TH1 *AddOverflowX(TH1 *h,double widthX) {
1145 // add overflow bin to X-axis
1146 int nx=h->GetNbinsX();
1147 double *xBins=new double[nx+2];
1148 for(int i=1;i<=nx;i++) {
1149 xBins[i-1]=h->GetXaxis()->GetBinLowEdge(i);
1150 }
1151 xBins[nx]=h->GetXaxis()->GetBinUpEdge(nx);
1152 xBins[nx+1]=xBins[nx]+widthX;
1153 TString name(h->GetName());
1154 name+="U";
1155 TH1 *r=new TH1D(name,h->GetTitle(),nx+1,xBins);
1156 for(int ix=0;ix<=nx+1;ix++) {
1157 r->SetBinContent(ix,h->GetBinContent(ix));
1158 r->SetBinError(ix,h->GetBinError(ix));
1159 }
1160 delete [] xBins;
1161 return r;
1162}
1163
1164void DrawOverflowX(TH1 *h,double posy) {
1165 double x1=h->GetXaxis()->GetBinLowEdge(h->GetNbinsX());
1166 double x2=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX());
1167 double y0=h->GetYaxis()->GetBinLowEdge(1);
1168 double y2=h->GetYaxis()->GetBinUpEdge(h->GetNbinsY());;
1169 if(h->GetDimension()==1) {
1170 y0=h->GetMinimum();
1171 y2=h->GetMaximum();
1172 }
1173 double w1=-0.3;
1174 TText *textX=new TText((1.+w1)*x2-w1*x1,(1.-posy)*y0+posy*y2,"Overflow bin");
1175 textX->SetNDC(kFALSE);
1176 textX->SetTextSize(0.05);
1177 textX->SetTextAngle(90.);
1178 textX->Draw();
1179 TLine *lineX=new TLine(x1,y0,x1,y2);
1180 lineX->Draw();
1181}
1182
1183void DrawOverflowY(TH1 *h,double posx) {
1184 double x0=h->GetXaxis()->GetBinLowEdge(1);
1185 double x2=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX());
1186 double y1=h->GetYaxis()->GetBinLowEdge(h->GetNbinsY());;
1187 double y2=h->GetYaxis()->GetBinUpEdge(h->GetNbinsY());;
1188 double w1=-0.3;
1189 TText *textY=new TText((1.-posx)*x0+posx*x2,(1.+w1)*y1-w1*y2,"Overflow bin");
1190 textY->SetNDC(kFALSE);
1191 textY->SetTextSize(0.05);
1192 textY->Draw();
1193 TLine *lineY=new TLine(x0,y1,x2,y1);
1194 lineY->Draw();
1195}
1196
1197void DrawPadProbability(TH2 *h) {
1198 h->Draw("COLZ");
1199 h->SetTitle("migration probabilities;P_{T}(gen) [GeV];P_{T}(rec) [GeV]");
1200 DrawOverflowX(h,0.05);
1201 DrawOverflowY(h,0.35);
1202}
1203
1204void DrawPadEfficiency(TH1 *h) {
1205 h->SetTitle("efficiency;P_{T}(gen) [GeV];#epsilon");
1206 h->SetLineColor(kBlue);
1207 h->SetMinimum(0.75);
1208 h->SetMaximum(1.0);
1209 h->Draw();
1210 DrawOverflowX(h,0.05);
1211 TLegend *legEfficiency=new TLegend(0.3,0.58,0.6,0.75);
1212 legEfficiency->SetBorderSize(0);
1213 legEfficiency->SetFillStyle(0);
1214 legEfficiency->SetTextSize(kLegendFontSize);
1215 legEfficiency->AddEntry(h,"reconstruction","l");
1216 legEfficiency->AddEntry((TObject*)nullptr," efficiency","");
1217 legEfficiency->Draw();
1218}
1219
1222 //gPad->SetLogy(kTRUE);
1223 double amax=0.0;
1224 for(int i=1;i<=histMcRec->GetNbinsX();i++) {
1225 amax=TMath::Max(amax,histMcRec->GetBinContent(i)
1226 +5.0*histMcRec->GetBinError(i));
1227 amax=TMath::Max(amax,histDataRec->GetBinContent(i)
1228 +2.0*histDataRec->GetBinError(i));
1229 }
1230 histMcRec->SetTitle("Reconstructed;P_{T}(rec);Nevent / GeV");
1231 histMcRec->Draw("HIST");
1232 histMcRec->SetLineColor(kBlue);
1233 histMcRec->SetMinimum(1.0);
1234 histMcRec->SetMaximum(amax);
1235 //histMcbgrRec->SetFillMode(1);
1236 histMcbgrRec->SetLineColor(kBlue-6);
1237 histMcbgrRec->SetFillColor(kBlue-10);
1238 histMcbgrRec->Draw("SAME HIST");
1239
1240 TH1 * histFoldBack=nullptr;
1242 histFoldBack=(TH1 *)
1243 histMcRec->Clone(histDataUnfold->GetName()+TString("_folded"));
1244 int nrec=histFoldBack->GetNbinsX();
1245 if((nrec==histProbability->GetNbinsY())&&
1246 (nrec==histMcbgrRec->GetNbinsX())&&
1247 (nrec==histDataRec->GetNbinsX())
1248 ) {
1249 for(int ix=1;ix<=nrec;ix++) {
1250 double sum=0.0;
1251 double sume2=0.0;
1252 for(int iy=0;iy<=histProbability->GetNbinsX()+1;iy++) {
1253 sum += histDataUnfold->GetBinContent(iy)*
1254 histDataUnfold->GetBinWidth(iy)*
1255 histProbability->GetBinContent(iy,ix);
1256 for(int iy2=0;iy2<=histProbability->GetNbinsX()+1;iy2++) {
1257 sume2 += histDataUnfold->GetBinError(iy)*
1258 histDataUnfold->GetBinWidth(iy)*
1259 histProbability->GetBinContent(iy,ix)*
1260 histDataUnfold->GetBinError(iy2)*
1261 histDataUnfold->GetBinWidth(iy2)*
1262 histProbability->GetBinContent(iy2,ix)*
1263 histRhoij->GetBinContent(iy,iy2);
1264 }
1265 }
1266 sum /= histFoldBack->GetBinWidth(ix);
1267 sum += histMcbgrRec->GetBinContent(ix);
1268 histFoldBack->SetBinContent(ix,sum);
1269 histFoldBack->SetBinError(ix,TMath::Sqrt(sume2)
1270 /histFoldBack->GetBinWidth(ix));
1271 }
1272 } else {
1273 cout<<"can not fold back: "<<nrec
1274 <<" "<<histProbability->GetNbinsY()
1275 <<" "<<histMcbgrRec->GetNbinsX()
1276 <<" "<<histDataRec->GetNbinsX()
1277 <<"\n";
1278 exit(0);
1279 }
1280
1281 histFoldBack->SetLineColor(kBlack);
1282 histFoldBack->SetMarkerStyle(0);
1283 histFoldBack->Draw("SAME HIST");
1284 }
1285
1286 histDataRec->SetLineColor(kRed);
1287 histDataRec->SetMarkerColor(kRed);
1288 histDataRec->Draw("SAME");
1290
1291 TLegend *legRec=new TLegend(0.4,0.5,0.68,0.85);
1292 legRec->SetBorderSize(0);
1293 legRec->SetFillStyle(0);
1294 legRec->SetTextSize(kLegendFontSize);
1295 legRec->AddEntry(histMcRec,"MC total","l");
1296 legRec->AddEntry(histMcbgrRec,"background","f");
1297 if(histFoldBack) {
1298 int ndf=-kNbinC;
1299 double sumD=0.,sumF=0.,chi2=0.;
1300 for(int i=1;i<=histDataRec->GetNbinsX();i++) {
1301 //cout<<histDataRec->GetBinContent(i)<<" "<<histFoldBack->GetBinContent(i)<<" "<<" w="<<histFoldBack->GetBinWidth(i)<<"\n";
1302 sumD+=histDataRec->GetBinContent(i)*histDataRec->GetBinWidth(i);
1303 sumF+=histFoldBack->GetBinContent(i)*histFoldBack->GetBinWidth(i);
1304 double pull=(histFoldBack->GetBinContent(i)-histDataRec->GetBinContent(i))/histDataRec->GetBinError(i);
1305 chi2+= pull*pull;
1306 ndf+=1;
1307 }
1308 legRec->AddEntry(histDataRec,TString::Format("data N_{evt}=%.0f",sumD),"lp");
1309 legRec->AddEntry(histFoldBack,TString::Format("folded N_{evt}=%.0f",sumF),"l");
1310 legRec->AddEntry((TObject*)nullptr,TString::Format("#chi^{2}=%.1f ndf=%d",chi2,ndf),"");
1311 //exit(0);
1312 } else {
1313 legRec->AddEntry(histDataRec,"data","lp");
1314 }
1315 legRec->Draw();
1316}
1317
1319 char const *text,double tau,vector<double> const *r,
1320 TF1 *f) {
1321 //gPad->SetLogy(kTRUE);
1322 double amin=0.;
1323 double amax=0.;
1324 for(int i=1;i<=histMcsigGen->GetNbinsX();i++) {
1325 if(histDataUnfold) {
1326 amin=TMath::Min(amin,histDataUnfold->GetBinContent(i)
1327 -1.1*histDataUnfold->GetBinError(i));
1328 amax=TMath::Max(amax,histDataUnfold->GetBinContent(i)
1329 +1.1*histDataUnfold->GetBinError(i));
1330 }
1331 amin=TMath::Min(amin,histMcsigGen->GetBinContent(i)
1332 -histMcsigGen->GetBinError(i));
1333 amin=TMath::Min(amin,histDataGen->GetBinContent(i)
1334 -histDataGen->GetBinError(i));
1335 amax=TMath::Max(amax,histMcsigGen->GetBinContent(i)
1336 +10.*histMcsigGen->GetBinError(i));
1337 amax=TMath::Max(amax,histDataGen->GetBinContent(i)
1338 +2.*histDataGen->GetBinError(i));
1339 }
1340 histMcsigGen->SetMinimum(amin);
1341 histMcsigGen->SetMaximum(amax);
1342
1343 histMcsigGen->SetTitle("Truth;P_{T};Nevent / GeV");
1344 histMcsigGen->SetLineColor(kBlue);
1345 histMcsigGen->Draw("HIST");
1346 histDataGen->SetLineColor(kRed);
1347 histDataGen->SetMarkerColor(kRed);
1348 histDataGen->SetMarkerSize(1.0);
1349 histDataGen->Draw("SAME HIST");
1350 if(histDataUnfold) {
1351 histDataUnfold->SetMarkerStyle(21);
1352 histDataUnfold->SetMarkerSize(0.7);
1353 histDataUnfold->Draw("SAME");
1354 }
1356
1357 if(f) {
1358 f->SetLineStyle(kSolid);
1359 f->Draw("SAME");
1360 }
1361
1362 TLegend *legTruth=new TLegend(0.32,0.65,0.6,0.9);
1363 legTruth->SetBorderSize(0);
1364 legTruth->SetFillStyle(0);
1365 legTruth->SetTextSize(kLegendFontSize);
1366 legTruth->AddEntry(histMcsigGen,"MC","l");
1367 if(!histDataUnfold) legTruth->AddEntry((TObject *)nullptr," Landau(5,2)","");
1368 legTruth->AddEntry(histDataGen,"data","l");
1369 if(!histDataUnfold) legTruth->AddEntry((TObject *)nullptr," Landau(6,1.8)","");
1370 if(histDataUnfold) {
1371 TString t;
1372 if(text) t=text;
1373 else t=histDataUnfold->GetName();
1374 if(tau>0) {
1375 t+=TString::Format(" #tau=%.2g",tau);
1376 }
1377 legTruth->AddEntry(histDataUnfold,t,"lp");
1378 if(r) {
1379 legTruth->AddEntry((TObject *)nullptr,"test wrt data:","");
1380 legTruth->AddEntry((TObject *)nullptr,TString::Format
1381 ("#chi^{2}/%d=%.1f prob=%.3f",
1382 (int)(*r)[3],(*r)[2]/(*r)[3],
1383 TMath::Prob((*r)[2],(*r)[3])),"");
1384 }
1385 }
1386 if(f) {
1387 legTruth->AddEntry(f,"fit","l");
1388 }
1389 legTruth->Draw();
1390 if(histDataUnfold ) {
1391 TPad *subpad = new TPad("subpad","",0.35,0.29,0.88,0.68);
1392 subpad->SetFillStyle(0);
1393 subpad->Draw();
1394 subpad->cd();
1395 amin=0.;
1396 amax=0.;
1397 int istart=11;
1398 for(int i=istart;i<=histMcsigGen->GetNbinsX();i++) {
1399 amin=TMath::Min(amin,histMcsigGen->GetBinContent(i)
1400 -histMcsigGen->GetBinError(i));
1401 amin=TMath::Min(amin,histDataGen->GetBinContent(i)
1402 -histDataGen->GetBinError(i));
1403 amin=TMath::Min(amin,histDataUnfold->GetBinContent(i)
1404 -histDataUnfold->GetBinError(i));
1405 amax=TMath::Max(amax,histMcsigGen->GetBinContent(i)
1406 +histMcsigGen->GetBinError(i));
1407 amax=TMath::Max(amax,histDataGen->GetBinContent(i)
1408 +histDataGen->GetBinError(i));
1409 amax=TMath::Max(amax,histDataUnfold->GetBinContent(i)
1410 +histDataUnfold->GetBinError(i));
1411 }
1412 TH1 *copyMcsigGen=(TH1*)histMcsigGen->Clone();
1413 TH1 *copyDataGen=(TH1*)histDataGen->Clone();
1415 copyMcsigGen->GetXaxis()->SetRangeUser
1416 (copyMcsigGen->GetXaxis()->GetBinLowEdge(istart),
1417 copyMcsigGen->GetXaxis()->GetBinUpEdge(copyMcsigGen->GetNbinsX()-1));
1418 copyMcsigGen->SetTitle(";;");
1419 copyMcsigGen->GetYaxis()->SetRangeUser(amin,amax);
1420 copyMcsigGen->Draw("HIST");
1421 copyDataGen->Draw("SAME HIST");
1422 copyDataUnfold->Draw("SAME");
1423 if(f) {
1424 ((TF1 *)f->Clone())->Draw("SAME");
1425 }
1426 }
1427}
1428
1430 vector<pair<TF1*,vector<double> > > const *table) {
1431 h->SetMinimum(-1.);
1432 h->SetMaximum(1.);
1433 h->SetTitle("correlation coefficients;P_{T}(gen) [GeV];P_{T}(gen) [GeV]");
1434 h->Draw("COLZ");
1435 DrawOverflowX(h,0.05);
1436 DrawOverflowY(h,0.05);
1437 if(table) {
1438 TLegend *legGCor=new TLegend(0.13,0.6,0.5,0.8);
1439 legGCor->SetBorderSize(0);
1440 legGCor->SetFillStyle(0);
1441 legGCor->SetTextSize(kLegendFontSize);
1442 vector<double> const &r=(*table)[table->size()-1].second;
1443 legGCor->AddEntry((TObject *)nullptr,TString::Format("min(#rho_{ij})=%5.2f",r[6]),"");
1444 legGCor->AddEntry((TObject *)nullptr,TString::Format("max(#rho_{ij})=%5.2f",r[7]),"");
1445 legGCor->AddEntry((TObject *)nullptr,TString::Format("avg(#rho_i)=%5.2f",r[5]),"");
1446 legGCor->Draw();
1447 }
1448}
1449
1450TH1 *g_fcnHist=nullptr;
1451TMatrixD *g_fcnMatrix=nullptr;
1452
1454 if(flag==0) {
1455 cout<<"fcn flag=0: npar="<<npar<<" gin="<<gin<<" par=[";
1456 for(int i=0;i<npar;i++) {
1457 cout<<" "<<u[i];
1458 }
1459 cout<<"]\n";
1460 }
1461 int n=g_fcnMatrix->GetNrows();
1462 TVectorD dy(n);
1463 double x0=0,y0=0.;
1464 for(int i=0;i<=n;i++) {
1465 double x1;
1466 if(i<1) x1=g_fcnHist->GetXaxis()->GetBinLowEdge(i+1);
1467 else x1=g_fcnHist->GetXaxis()->GetBinUpEdge(i);
1468 double y1=TMath::LandauI((x1-u[1])/u[2]);
1469 if(i>0) {
1470 double iy=u[0]*u[2]*(y1-y0)/(x1-x0);
1471 dy(i-1)=iy-g_fcnHist->GetBinContent(i);
1472 //cout<<"i="<<i<<" iy="<<iy<<" delta="<< dy(i-1)<<"\n";
1473 }
1474 x0=x1;
1475 y0=y1;
1476 //cout<<"i="<<i<<" y1="<<y1<<" x1="<<x1<<"\n";
1477 }
1478 TVectorD Hdy=(*g_fcnMatrix) * dy;
1479 //Hdy.Print();
1480 f=Hdy*dy;
1481 //exit(0);
1482}
1483
1484
1485TFitResultPtr DoFit(TH1 *h,TH2 *rho,TH1 *truth,const char *text,
1486 vector<pair<TF1 *,vector<double> > > &table,int niter) {
1487 TString option="IESN";
1488 cout<<h->GetName()<<"\n";
1489 double gcorAvg=0.;
1490 double rhoMin=0.;
1491 double rhoMax=0.;
1492 if(rho) {
1493 g_fcnHist=h;
1494 int n=h->GetNbinsX()-1; // overflow is included as extra bin, exclude in fit
1495 TMatrixDSym v(n);
1496 //g_fcnMatrix=new TMatrixD(n,n);
1497 for(int i=0;i<n;i++) {
1498 for(int j=0;j<n;j++) {
1499 v(i,j)=rho->GetBinContent(i+1,j+1)*
1500 (h->GetBinError(i+1)*h->GetBinError(j+1));
1501 }
1502 }
1504 TMatrixD d(n,n);
1505 TVectorD di(ev.GetEigenValues());
1506 for(int i=0;i<n;i++) {
1507 if(di(i)>0.0) {
1508 d(i,i)=1./di(i);
1509 } else {
1510 cout<<"bad eigenvalue i="<<i<<" di="<<di(i)<<"\n";
1511 exit(0);
1512 }
1513 }
1514 TMatrixD O(ev.GetEigenVectors());
1518 int error=0;
1519 for(int i=0;i<n;i++) {
1520 if(TMath::Abs(test(i,i)-1.0)>1.E-7) {
1521 error++;
1522 }
1523 for(int j=0;j<n;j++) {
1524 if(i==j) continue;
1525 if(TMath::Abs(test(i,j)>1.E-7)) error++;
1526 }
1527 }
1528 // calculate global correlation coefficient (all bins)
1529 TMatrixDSym v1(n+1);
1530 rhoMin=1.;
1531 rhoMax=-1.;
1532 for(int i=0;i<=n;i++) {
1533 for(int j=0;j<=n;j++) {
1534 double rho_ij=rho->GetBinContent(i+1,j+1);
1535 v1(i,j)=rho_ij*
1536 (h->GetBinError(i+1)*h->GetBinError(j+1));
1537 if(i!=j) {
1540 }
1541 }
1542 }
1544 TMatrixD d1(n+1,n+1);
1545 TVectorD di1(ev1.GetEigenValues());
1546 for(int i=0;i<=n;i++) {
1547 if(di1(i)>0.0) {
1548 d1(i,i)=1./di1(i);
1549 } else {
1550 cout<<"bad eigenvalue i="<<i<<" di1="<<di1(i)<<"\n";
1551 exit(0);
1552 }
1553 }
1554 TMatrixD O1(ev1.GetEigenVectors());
1557 for(int i=0;i<=n;i++) {
1558 double gcor2=1.-1./(vinv1(i,i)*v1(i,i));
1559 if(gcor2>=0.0) {
1560 double gcor=TMath::Sqrt(gcor2);
1561 gcorAvg += gcor;
1562 } else {
1563 cout<<"bad global correlation "<<i<<" "<<gcor2<<"\n";
1564 }
1565 }
1566 gcorAvg /=(n+1);
1567 /* if(error) {
1568 v.Print();
1569 g_fcnMatrix->Print();
1570 exit(0);
1571 } */
1572 //g_fcnMatrix->Invert();
1573 //from: HFitImpl.cxx
1574 // TVirtualFitter::FCNFunc_t userFcn = 0;
1575 // typedef void (* FCNFunc_t )(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
1576 // userFcn = (TVirtualFitter::GetFitter())->GetFCN();
1577 // (TVirtualFitter::GetFitter())->SetUserFunc(f1);
1578 //...
1579 //fitok = fitter->FitFCN( userFcn );
1580
1581 TVirtualFitter::Fitter(h)->SetFCN(fcn);
1582 option += "U";
1583
1584 }
1585 double xmax=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX()-1);
1586 TF1 *landau=new TF1(text,"[0]*TMath::Landau(x,[1],[2],0)",
1587 0.,xmax);
1588 landau->SetParameter(0,6000.);
1589 landau->SetParameter(1,5.);
1590 landau->SetParameter(2,2.);
1591 landau->SetParError(0,10.);
1592 landau->SetParError(1,0.5);
1593 landau->SetParError(2,0.1);
1594 TFitResultPtr s=h->Fit(landau,option,nullptr,0.,xmax);
1595 vector<double> r(8);
1596 int np=landau->GetNpar();
1597 fcn(np,nullptr,r[0],landau->GetParameters(),0);
1598 r[1]=h->GetNbinsX()-1-landau->GetNpar();
1599 for(int i=0;i<h->GetNbinsX()-1;i++) {
1600 double di=h->GetBinContent(i+1)-truth->GetBinContent(i+1);
1601 if(g_fcnMatrix) {
1602 for(int j=0;j<h->GetNbinsX()-1;j++) {
1603 double dj=h->GetBinContent(j+1)-truth->GetBinContent(j+1);
1604 r[2]+=di*dj*(*g_fcnMatrix)(i,j);
1605 }
1606 } else {
1607 double pull=di/h->GetBinError(i+1);
1608 r[2]+=pull*pull;
1609 }
1610 r[3]+=1.0;
1611 }
1612 r[4]=niter;
1613 if(!niter) r[4]=0.25;
1614 r[5]=gcorAvg;
1615 r[6]=rhoMin;
1616 r[7]=rhoMax;
1617 if(rho) {
1618 g_fcnHist=nullptr;
1619 delete g_fcnMatrix;
1620 g_fcnMatrix=nullptr;
1621 }
1622 table.push_back(make_pair(landau,r));
1623 return s;
1624}
1625
1626#ifdef WITH_IDS
1627
1628//===================== interface to IDS unfolding code follows here
1629// contact Bogdan Malescu to find it
1630
1631#include "ids_code.cc"
1632
1634
1635 int N_=data->GetNrows();
1636 soustr = new TVectorD(N_);
1637 for( Int_t i=0; i<N_; i++ ){ (*soustr)[i] = 0.; }
1638 unfres1IDS_ = Unfold( data, dataErr, A_, N_, lambdaL_, soustr );
1639}
1640
1642
1643 int N_=data->GetNrows();
1645 delete unfres2IDS_;
1646 unfres2IDS_ = Unfold( data, dataErr, Am_, N_, lambdaU_, soustr );
1647}
1648
1649#endif
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
@ kRed
Definition Rtypes.h:67
@ kOrange
Definition Rtypes.h:68
@ kBlack
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:67
@ kWhite
Definition Rtypes.h:66
@ kCyan
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
@ kAzure
Definition Rtypes.h:68
@ kSolid
Definition TAttLine.h:52
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
constexpr Int_t kInfo
Definition TError.h:45
Int_t gErrorIgnoreLevel
errors with level below this value will be ignored. Default is kUnset.
Definition TError.cxx:33
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t width
Option_t Option_t style
Option_t Option_t TPoint TPoint const char text
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
Double_t GetBinError(Int_t bin) const override
float xmax
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
#define gPad
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:40
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:182
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
A TGraphErrors is a TGraph with error bars.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:927
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
Definition TH1.cxx:9244
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
When this static function is called with sumw2=kTRUE, all new histograms will automatically activate ...
Definition TH1.cxx:6751
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9260
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:400
Service class for 2-D histogram classes.
Definition TH2.h:30
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:97
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
Use the TLine constructor to create a simple line.
Definition TLine.h:22
TMatrixDSymEigen.
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
The most important graphics class in the ROOT system.
Definition TPad.h:28
Random number generator class based on M.
Definition TRandom3.h:27
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Basic string class.
Definition TString.h:138
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1641
Base class for several text objects.
Definition TText.h:22
Binning schemes for use with the unfolding algorithm TUnfoldDensity.
TH1 * CreateHistogram(const char *histogramName, Bool_t originalAxisBinning=kFALSE, Int_t **binMap=nullptr, const char *histogramTitle=nullptr, const char *axisSteering=nullptr) const
create a THxx histogram capable to hold the bins of this binning node and its children
Double_t GetBinSize(Int_t iBin) const
get N-dimensional bin size
Int_t GetEndBin(void) const
last+1 bin of this node (includes children)
An algorithm to unfold distributions from detector to truth level.
@ kEScanTauRhoAvg
average global correlation coefficient (from TUnfold::GetRhoI())
@ kDensityModeNone
no scale factors, matrix L is similar to unity matrix
@ kEConstraintArea
enforce preservation of the area
Definition TUnfold.h:119
@ kEConstraintNone
use no extra constraint
Definition TUnfold.h:116
ERegMode
choice of regularisation scheme
Definition TUnfold.h:123
@ kRegModeSize
regularise the amplitude of the output distribution
Definition TUnfold.h:129
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition TUnfold.h:146
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
TLine * line
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
double landau(double x, double mu, double sigma)
Definition MathFuncs.h:374
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251
Double_t Prob(Double_t chi2, Int_t ndf)
Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf...
Definition TMath.cxx:637
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:732
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:199
Double_t Hypot(Double_t x, Double_t y)
Returns sqrt(x*x + y*y)
Definition TMath.cxx:59
Double_t LandauI(Double_t x)
Returns the cumulative (lower tail integral) of the Landau distribution function at point x.
Definition TMath.cxx:2847
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
th1 Draw()
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2339