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