Logo ROOT   6.07/09
Reference Guide
rf602_chi2fit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #602
5 ///
6 /// Setting up a chi^2 fit to a binned dataset
7 ///
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "RooConstVar.h"
17 #include "RooChebychev.h"
18 #include "RooAddPdf.h"
19 #include "RooChi2Var.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 using namespace RooFit ;
24 
25 
26 void rf602_chi2fit()
27 {
28 
29  // S e t u p m o d e l
30  // ---------------------
31 
32  // Declare observable x
33  RooRealVar x("x","x",0,10) ;
34 
35  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
36  RooRealVar mean("mean","mean of gaussians",5) ;
37  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
38  RooRealVar sigma2("sigma2","width of gaussians",1) ;
39 
40  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
41  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
42 
43  // Build Chebychev polynomial p.d.f.
44  RooRealVar a0("a0","a0",0.5,0.,1.) ;
45  RooRealVar a1("a1","a1",0.2,0.,1.) ;
46  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
47 
48  // Sum the signal components into a composite signal p.d.f.
49  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
50  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
51 
52  // Sum the composite signal and background
53  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
54  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
55 
56 
57  // C r e a t e b i n n e d d a t a s e t
58  // -----------------------------------------
59 
60  RooDataSet* d = model.generate(x,10000) ;
61  RooDataHist* dh = d->binnedClone() ;
62 
63  // Construct a chi^2 of the data and the model.
64  // When a p.d.f. is used in a chi^2 fit, the probability density scaled
65  // by the number of events in the dataset to obtain the fit function
66  // If model is an extended p.d.f, the expected number events is used
67  // instead of the observed number of events.
68  model.chi2FitTo(*dh) ;
69 
70  // NB: It is also possible to fit a RooAbsReal function to a RooDataHist
71  // using chi2FitTo().
72 
73  // Note that entries with zero bins are _not_ allowed
74  // for a proper chi^2 calculation and will give error
75  // messages
76  RooDataSet* dsmall = (RooDataSet*) d->reduce(EventRange(1,100)) ;
77  RooDataHist* dhsmall = dsmall->binnedClone() ;
78  RooChi2Var chi2_lowstat("chi2_lowstat","chi2",model,*dhsmall) ;
79  cout << chi2_lowstat.getVal() << endl ;
80 
81 
82 }
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooAbsData * reduce(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg())
Create a reduced copy of this dataset.
Definition: RooAbsData.cxx:343
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooCmdArg EventRange(Int_t nStart, Int_t nStop)
RooDataHist * binnedClone(const char *newName=0, const char *newTitle=0) const
Return binned clone of this dataset.
Definition: RooDataSet.cxx:981
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25