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