Logo ROOT  
Reference Guide
rf313_paramranges.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## Multidimensional models: working with parameterized ranges to define non-rectangular regions for fitting and integration
6##
7## \macro_code
8##
9## \date February 2018
10## \author Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14
15# Create 3D pdf
16# -------------------------
17
18# Define observable (x,y,z)
19x = ROOT.RooRealVar("x", "x", 0, 10)
20y = ROOT.RooRealVar("y", "y", 0, 10)
21z = ROOT.RooRealVar("z", "z", 0, 10)
22
23# Define 3 dimensional pdf
24z0 = ROOT.RooRealVar("z0", "z0", -0.1, 1)
25px = ROOT.RooPolynomial(
26 "px", "px", x, ROOT.RooArgList(
27 ROOT.RooFit.RooConst(0)))
28py = ROOT.RooPolynomial(
29 "py", "py", y, ROOT.RooArgList(
30 ROOT.RooFit.RooConst(0)))
31pz = ROOT.RooPolynomial("pz", "pz", z, ROOT.RooArgList(z0))
32pxyz = ROOT.RooProdPdf("pxyz", "pxyz", ROOT.RooArgList(px, py, pz))
33
34# Defined non-rectangular region R in (x, y, z)
35# -------------------------------------------------------------------------------------
36
37#
38# R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
39#
40
41# Construct range parameterized in "R" in y [ 0.1*x, 0.9*x ]
42ylo = ROOT.RooFormulaVar("ylo", "0.1*x", ROOT.RooArgList(x))
43yhi = ROOT.RooFormulaVar("yhi", "0.9*x", ROOT.RooArgList(x))
44y.setRange("R", ylo, yhi)
45
46# Construct parameterized ranged "R" in z [ 0, 0.1*y^2 ]
47zlo = ROOT.RooFormulaVar("zlo", "0.0*y", ROOT.RooArgList(y))
48zhi = ROOT.RooFormulaVar("zhi", "0.1*y*y", ROOT.RooArgList(y))
49z.setRange("R", zlo, zhi)
50
51# Calculate integral of normalized pdf in R
52# ----------------------------------------------------------------------------------
53
54# Create integral over normalized pdf model over x,y, in "R" region
55intPdf = pxyz.createIntegral(ROOT.RooArgSet(
56 x, y, z), ROOT.RooArgSet(x, y, z), "R")
57
58# Plot value of integral as function of pdf parameter z0
59frame = z0.frame(ROOT.RooFit.Title(
60 "Integral of pxyz over x,y, in region R"))
61intPdf.plotOn(frame)
62
63c = ROOT.TCanvas("rf313_paramranges", "rf313_paramranges", 600, 600)
64ROOT.gPad.SetLeftMargin(0.15)
65frame.GetYaxis().SetTitleOffset(1.6)
66frame.Draw()
67
68c.SaveAs("rf313_paramranges.png")