ROOT
master
Reference Guide
Loading...
Searching...
No Matches
rf313_paramranges.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
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
14
import
ROOT
15
16
17
# Create 3D pdf
18
# -------------------------
19
20
# Define observable (x,y,z)
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, 0, 10)
22
y =
ROOT.RooRealVar
(
"y"
,
"y"
, 0, 10)
23
z =
ROOT.RooRealVar
(
"z"
,
"z"
, 0, 10)
24
25
# Define 3 dimensional pdf
26
z0 =
ROOT.RooRealVar
(
"z0"
,
"z0"
, -0.1, 1)
27
px =
ROOT.RooPolynomial
(
"px"
,
"px"
, x, [0.0])
28
py =
ROOT.RooPolynomial
(
"py"
,
"py"
, y, [0.0])
29
pz =
ROOT.RooPolynomial
(
"pz"
,
"pz"
, z, [z0])
30
pxyz =
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 ]
40
ylo =
ROOT.RooFormulaVar
(
"ylo"
,
"0.1*x"
, [x])
41
yhi =
ROOT.RooFormulaVar
(
"yhi"
,
"0.9*x"
, [x])
42
y.setRange
(
"R"
, ylo, yhi)
43
44
# Construct parameterized ranged "R" in z [ 0, 0.1*y^2 ]
45
zlo =
ROOT.RooFormulaVar
(
"zlo"
,
"0.0*y"
, [y])
46
zhi =
ROOT.RooFormulaVar
(
"zhi"
,
"0.1*y*y"
, [y])
47
z.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
53
intPdf =
pxyz.createIntegral
({x, y, z}, {x, y, z},
"R"
)
54
55
# Plot value of integral as function of pdf parameter z0
56
frame =
z0.frame
(Title=
"Integral of pxyz over x,y, in region R"
)
57
intPdf.plotOn
(frame)
58
59
c =
ROOT.TCanvas
(
"rf313_paramranges"
,
"rf313_paramranges"
, 600, 600)
60
ROOT.gPad.SetLeftMargin
(0.15)
61
frame.GetYaxis
().SetTitleOffset(1.6)
62
frame.Draw
()
63
64
c.SaveAs
(
"rf313_paramranges.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
roofit
rf313_paramranges.py
ROOT master - Reference Guide Generated on Fri Jan 24 2025 10:31:07 (GVA Time) using Doxygen 1.10.0