Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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(
23 "gx", "gx", x, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
24gy = ROOT.RooGaussian(
25 "gy", "gy", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
26gz = ROOT.RooGaussian(
27 "gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
28sig = ROOT.RooProdPdf("sig", "sig", ROOT.RooArgList(gx, gy, gz))
29
30# Create background pdf poly(x)*poly(y)*poly(z)
31px = ROOT.RooPolynomial("px", "px", x, ROOT.RooArgList(
32 ROOT.RooFit.RooConst(-0.1), ROOT.RooFit.RooConst(0.004)))
33py = ROOT.RooPolynomial("py", "py", y, ROOT.RooArgList(
34 ROOT.RooFit.RooConst(0.1), ROOT.RooFit.RooConst(-0.004)))
35pz = ROOT.RooPolynomial("pz", "pz", z)
36bkg = ROOT.RooProdPdf("bkg", "bkg", ROOT.RooArgList(px, py, pz))
37
38# Create composite pdf sig+bkg
39fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0., 1.)
40model = ROOT.RooAddPdf(
41 "model", "model", ROOT.RooArgList(
42 sig, bkg), ROOT.RooArgList(fsig))
43
44data = model.generate(ROOT.RooArgSet(x, y, z), 20000)
45
46# Project pdf and data on x
47# -------------------------------------------------
48
49# Make plain projection of data and pdf on x observable
50frame = x.frame(ROOT.RooFit.Title(
51 "Projection of 3D data and pdf on X"), ROOT.RooFit.Bins(40))
52data.plotOn(frame)
53model.plotOn(frame)
54
55# Project pdf and data on x in signal range
56# ----------------------------------------------------------------------------------
57
58# Define signal region in y and z observables
59y.setRange("sigRegion", -1, 1)
60z.setRange("sigRegion", -1, 1)
61
62# Make plot frame
63frame2 = x.frame(ROOT.RooFit.Title(
64 "Same projection on X in signal range of (Y,Z)"), ROOT.RooFit.Bins(40))
65
66# Plot subset of data in which all observables are inside "sigRegion"
67# For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
68# an implicit definition is used that is identical to the full range (i.e.
69# [-5,5] for x)
70data.plotOn(frame2, ROOT.RooFit.CutRange("sigRegion"))
71
72# Project model on x, projected observables (y,z) only in "sigRegion"
73model.plotOn(frame2, ROOT.RooFit.ProjectionRange("sigRegion"))
74
75c = ROOT.TCanvas("rf311_rangeplot", "rf310_rangeplot", 800, 400)
76c.Divide(2)
77c.cd(1)
78ROOT.gPad.SetLeftMargin(0.15)
79frame.GetYaxis().SetTitleOffset(1.4)
80frame.Draw()
81c.cd(2)
82ROOT.gPad.SetLeftMargin(0.15)
83frame2.GetYaxis().SetTitleOffset(1.4)
84frame2.Draw()
85
86c.SaveAs("rf311_rangeplot.png")