Logo ROOT  
Reference Guide
rf311_rangeplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## Multidimensional models: projecting p.d.f and data ranges in continuous observables
6##
7## \macro_code
8##
9## \date February 2018
10## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14# Create 3D pdf and data
15# -------------------------------------------
16
17# Create observables
18x = ROOT.RooRealVar("x", "x", -5, 5)
19y = ROOT.RooRealVar("y", "y", -5, 5)
20z = ROOT.RooRealVar("z", "z", -5, 5)
21
22# Create signal pdf gauss(x)*gauss(y)*gauss(z)
23gx = ROOT.RooGaussian(
24 "gx", "gx", x, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
25gy = ROOT.RooGaussian(
26 "gy", "gy", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
27gz = ROOT.RooGaussian(
28 "gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
29sig = ROOT.RooProdPdf("sig", "sig", ROOT.RooArgList(gx, gy, gz))
30
31# Create background pdf poly(x)*poly(y)*poly(z)
32px = ROOT.RooPolynomial("px", "px", x, ROOT.RooArgList(
33 ROOT.RooFit.RooConst(-0.1), ROOT.RooFit.RooConst(0.004)))
34py = ROOT.RooPolynomial("py", "py", y, ROOT.RooArgList(
35 ROOT.RooFit.RooConst(0.1), ROOT.RooFit.RooConst(-0.004)))
36pz = ROOT.RooPolynomial("pz", "pz", z)
37bkg = ROOT.RooProdPdf("bkg", "bkg", ROOT.RooArgList(px, py, pz))
38
39# Create composite pdf sig+bkg
40fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0., 1.)
41model = ROOT.RooAddPdf(
42 "model", "model", ROOT.RooArgList(
43 sig, bkg), ROOT.RooArgList(fsig))
44
45data = model.generate(ROOT.RooArgSet(x, y, z), 20000)
46
47# Project pdf and data on x
48# -------------------------------------------------
49
50# Make plain projection of data and pdf on x observable
51frame = x.frame(ROOT.RooFit.Title(
52 "Projection of 3D data and pdf on X"), ROOT.RooFit.Bins(40))
53data.plotOn(frame)
54model.plotOn(frame)
55
56# Project pdf and data on x in signal range
57# ----------------------------------------------------------------------------------
58
59# Define signal region in y and z observables
60y.setRange("sigRegion", -1, 1)
61z.setRange("sigRegion", -1, 1)
62
63# Make plot frame
64frame2 = x.frame(ROOT.RooFit.Title(
65 "Same projection on X in signal range of (Y,Z)"), ROOT.RooFit.Bins(40))
66
67# Plot subset of data in which all observables are inside "sigRegion"
68# For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
69# an implicit definition is used that is identical to the full range (i.e.
70# [-5,5] for x)
71data.plotOn(frame2, ROOT.RooFit.CutRange("sigRegion"))
72
73# Project model on x, projected observables (y,z) only in "sigRegion"
74model.plotOn(frame2, ROOT.RooFit.ProjectionRange("sigRegion"))
75
76c = ROOT.TCanvas("rf311_rangeplot", "rf310_rangeplot", 800, 400)
77c.Divide(2)
78c.cd(1)
79ROOT.gPad.SetLeftMargin(0.15)
80frame.GetYaxis().SetTitleOffset(1.4)
81frame.Draw()
82c.cd(2)
83ROOT.gPad.SetLeftMargin(0.15)
84frame2.GetYaxis().SetTitleOffset(1.4)
85frame2.Draw()
86
87c.SaveAs("rf311_rangeplot.png")