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