Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf311_rangeplot.py File Reference

Namespaces

namespace  rf311_rangeplot
 

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: projecting pdf and data ranges in continuous observables

import ROOT
# Create 3D pdf and data
# -------------------------------------------
# Create observables
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -5, 5)
z = ROOT.RooRealVar("z", "z", -5, 5)
# Create signal pdf gauss(x)*gauss(y)*gauss(z)
gx = ROOT.RooGaussian("gx", "gx", x, 0.0, 1.0)
gy = ROOT.RooGaussian("gy", "gy", y, 0.0, 1.0)
gz = ROOT.RooGaussian("gz", "gz", z, 0.0, 1.0)
sig = ROOT.RooProdPdf("sig", "sig", [gx, gy, gz])
# Create background pdf poly(x)*poly(y)*poly(z)
px = ROOT.RooPolynomial("px", "px", x, [-0.1, 0.004])
py = ROOT.RooPolynomial("py", "py", y, [0.1, -0.004])
pz = ROOT.RooPolynomial("pz", "pz", z)
bkg = ROOT.RooProdPdf("bkg", "bkg", [px, py, pz])
# Create composite pdf sig+bkg
fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "model", [sig, bkg], [fsig])
data = model.generate({x, y, z}, 20000)
# Project pdf and data on x
# -------------------------------------------------
# Make plain projection of data and pdf on x observable
frame = x.frame(Title="Projection of 3D data and pdf on X", Bins=40)
data.plotOn(frame)
model.plotOn(frame)
# Project pdf and data on x in signal range
# ----------------------------------------------------------------------------------
# Define signal region in y and z observables
y.setRange("sigRegion", -1, 1)
z.setRange("sigRegion", -1, 1)
# Make plot frame
frame2 = x.frame(Title="Same projection on X in signal range of (Y,Z)", Bins=40)
# Plot subset of data in which all observables are inside "sigRegion"
# For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
# an implicit definition is used that is identical to the full range (i.e.
# [-5,5] for x)
data.plotOn(frame2, CutRange="sigRegion")
# Project model on x, projected observables (y,z) only in "sigRegion"
model.plotOn(frame2, ProjectionRange="sigRegion")
c = ROOT.TCanvas("rf311_rangeplot", "rf310_rangeplot", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.4)
frame2.Draw()
c.SaveAs("rf311_rangeplot.png")
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on x integrates over variables (z,y)
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'sigRegion' created with bounds [-1,1]
[#1] INFO:Eval -- RooRealVar::setRange(z) new range named 'sigRegion' created with bounds [-1,1]
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 1692 events out of 20000 total events
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on x integrates over variables (z,y) in range sigRegion
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf311_rangeplot.py.