Logo ROOT   6.16/01
Reference Guide
rf204_extrangefit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #204
5## Extended maximum likelihood fit with alternate range definition
6## for observed number of events.
7##
8## \macro_code
9##
10## \date February 2018
11## \author Clemens Lange
12## \author Wouter Verkerke (C version)
13
14import ROOT
15
16# Set up component pdfs
17# ---------------------------------------
18
19# Declare observable x
20x = ROOT.RooRealVar("x", "x", 0, 10)
21
22# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
23# their parameters
24mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
25sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
26sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
27
28sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
29sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
30
31# Build Chebychev polynomial p.d.f.
32a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
33a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
34bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
35
36# Sum the signal components into a composite signal p.d.f.
37sig1frac = ROOT.RooRealVar(
38 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
39sig = ROOT.RooAddPdf(
40 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
41
42# Construct extended comps with range spec
43# ------------------------------------------------------------------------------
44
45# Define signal range in which events counts are to be defined
46x.setRange("signalRange", 4, 6)
47
48# Associated nsig/nbkg as expected number of events with sig/bkg
49# _in_the_range_ "signalRange"
50nsig = ROOT.RooRealVar(
51 "nsig", "number of signal events in signalRange", 500, 0., 10000)
52nbkg = ROOT.RooRealVar(
53 "nbkg", "number of background events in signalRange", 500, 0, 10000)
54esig = ROOT.RooExtendPdf(
55 "esig", "extended signal p.d.f", sig, nsig, "signalRange")
56ebkg = ROOT.RooExtendPdf(
57 "ebkg", "extended background p.d.f", bkg, nbkg, "signalRange")
58
59# Sum extended components
60# ---------------------------------------------
61
62# Construct sum of two extended p.d.f. (no coefficients required)
63model = ROOT.RooAddPdf("model", "(g1+g2)+a", ROOT.RooArgList(ebkg, esig))
64
65# Sample data, fit model
66# -------------------------------------------
67
68# Generate 1000 events from model so that nsig, come out to numbers <<500
69# in fit
70data = model.generate(ROOT.RooArgSet(x), 1000)
71
72# Perform unbinned extended ML fit to data
73r = model.fitTo(data, ROOT.RooFit.Extended(ROOT.kTRUE), ROOT.RooFit.Save())
74r.Print()