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