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