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