Logo ROOT   6.16/01
Reference Guide
rf312_multirangefit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4##
5## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #312
6##
7## Performing fits in multiple (disjoint) ranges in one or more dimensions
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15import ROOT
16
17
18# Create 2D pdf and data
19# -------------------------------------------
20
21# Define observables x,y
22x = ROOT.RooRealVar("x", "x", -10, 10)
23y = ROOT.RooRealVar("y", "y", -10, 10)
24
25# Construct the signal pdf gauss(x)*gauss(y)
26mx = ROOT.RooRealVar("mx", "mx", 1, -10, 10)
27my = ROOT.RooRealVar("my", "my", 1, -10, 10)
28
29gx = ROOT.RooGaussian("gx", "gx", x, mx, ROOT.RooFit.RooConst(1))
30gy = ROOT.RooGaussian("gy", "gy", y, my, ROOT.RooFit.RooConst(1))
31
32sig = ROOT.RooProdPdf("sig", "sig", gx, gy)
33
34# Construct the background pdf (flat in x,y)
35px = ROOT.RooPolynomial("px", "px", x)
36py = ROOT.RooPolynomial("py", "py", y)
37bkg = ROOT.RooProdPdf("bkg", "bkg", px, py)
38
39# Construct the composite model sig+bkg
40f = ROOT.RooRealVar("f", "f", 0., 1.)
41model = ROOT.RooAddPdf("model", "model", ROOT.RooArgList(sig, bkg), ROOT.RooArgList(f))
42
43# Sample 10000 events in (x,y) from the model
44modelData = model.generate(ROOT.RooArgSet(x, y), 10000)
45
46# Define signal and sideband regions
47# -------------------------------------------------------------------
48
49# Construct the SideBand1,SideBand2, regions
50#
51# |
52# +-------------+-----------+
53# | | |
54# | Side | Sig |
55# | Band1 | nal |
56# | | |
57# --+-------------+-----------+--
58# | |
59# | Side |
60# | Band2 |
61# | |
62# +-------------+-----------+
63# |
64
65x.setRange("SB1", -10, +10)
66y.setRange("SB1", -10, 0)
67
68x.setRange("SB2", -10, 0)
69y.setRange("SB2", 0, +10)
70
71x.setRange("SIG", 0, +10)
72y.setRange("SIG", 0, +10)
73
74x.setRange("FULL", -10, +10)
75y.setRange("FULL", -10, +10)
76
77# Perform fits in individual sideband regions
78# -------------------------------------------------------------------------------------
79
80# Perform fit in SideBand1 region (ROOT.RooAddPdf coefficients will be
81# interpreted in full range)
82r_sb1 = model.fitTo(modelData, ROOT.RooFit.Range(
83 "SB1"), ROOT.RooFit.Save())
84
85# Perform fit in SideBand2 region (ROOT.RooAddPdf coefficients will be
86# interpreted in full range)
87r_sb2 = model.fitTo(modelData, ROOT.RooFit.Range(
88 "SB2"), ROOT.RooFit.Save())
89
90# Perform fits in joint sideband regions
91# -----------------------------------------------------------------------------
92
93# Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
94# (ROOT.RooAddPdf coefficients will be interpreted in full range)
95r_sb12 = model.fitTo(modelData, ROOT.RooFit.Range(
96 "SB1,SB2"), ROOT.RooFit.Save())
97
98# Print results for comparison
99r_sb1.Print()
100r_sb2.Print()
101r_sb12.Print()