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