Logo ROOT  
Reference Guide
rf602_chi2fit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #602
6##
7## Setting up a chi^2 fit to a binned dataset
8##
9## \macro_code
10##
11## \date February 2018
12## \authors Clemens Lange, Wouter Verkerke (C version)
13
14from __future__ import print_function
15import ROOT
16
17
18# Set up model
19# ---------------------
20
21# Declare observable x
22x = ROOT.RooRealVar("x", "x", 0, 10)
23
24# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
25# their parameters
26mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
27sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
28sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
29
30sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
31sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
32
33# Build Chebychev polynomial p.d.f.
34a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
35a1 = ROOT.RooRealVar("a1", "a1", 0.2, 0.0, 1.0)
36bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
37
38# Sum the signal components into a composite signal p.d.f.
39sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
40sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
41
42# Sum the composite signal and background
43bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
44model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
45
46# Create biuned dataset
47# -----------------------------------------
48
49d = model.generate({x}, 10000)
50dh = d.binnedClone()
51
52# Construct a chi^2 of the data and the model.
53# When a p.d.f. is used in a chi^2 fit, probability density scaled
54# by the number of events in the dataset to obtain the fit function
55# If model is an extended p.d.f, expected number events is used
56# instead of the observed number of events.
57ll = ROOT.RooLinkedList()
58model.chi2FitTo(dh, ll)
59
60# NB: It is also possible to fit a ROOT.RooAbsReal function to a ROOT.RooDataHist
61# using chi2FitTo().
62
63# Note that entries with zero bins are _not_ allowed
64# for a proper chi^2 calculation and will give error
65# messages
66dsmall = d.reduce(ROOT.RooFit.EventRange(1, 100))
67dhsmall = dsmall.binnedClone()
68chi2_lowstat = model.createChi2(dhsmall)
69print(chi2_lowstat.getVal())