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