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