ROOT
master
Reference Guide
Loading...
Searching...
No Matches
rf311_rangeplot.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
3
## \notebook
4
## Multidimensional models: projecting pdf and data ranges in continuous observables
5
##
6
## \macro_image
7
## \macro_code
8
## \macro_output
9
##
10
## \date February 2018
11
## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13
import
ROOT
14
15
# Create 3D pdf and data
16
# -------------------------------------------
17
18
# Create observables
19
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -5, 5)
20
y =
ROOT.RooRealVar
(
"y"
,
"y"
, -5, 5)
21
z =
ROOT.RooRealVar
(
"z"
,
"z"
, -5, 5)
22
23
# Create signal pdf gauss(x)*gauss(y)*gauss(z)
24
gx =
ROOT.RooGaussian
(
"gx"
,
"gx"
, x, 0.0, 1.0)
25
gy =
ROOT.RooGaussian
(
"gy"
,
"gy"
, y, 0.0, 1.0)
26
gz =
ROOT.RooGaussian
(
"gz"
,
"gz"
, z, 0.0, 1.0)
27
sig =
ROOT.RooProdPdf
(
"sig"
,
"sig"
, [gx, gy, gz])
28
29
# Create background pdf poly(x)*poly(y)*poly(z)
30
px =
ROOT.RooPolynomial
(
"px"
,
"px"
, x, [-0.1, 0.004])
31
py =
ROOT.RooPolynomial
(
"py"
,
"py"
, y, [0.1, -0.004])
32
pz =
ROOT.RooPolynomial
(
"pz"
,
"pz"
, z)
33
bkg =
ROOT.RooProdPdf
(
"bkg"
,
"bkg"
, [px, py, pz])
34
35
# Create composite pdf sig+bkg
36
fsig =
ROOT.RooRealVar
(
"fsig"
,
"signal fraction"
, 0.1, 0.0, 1.0)
37
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [sig, bkg], [fsig])
38
39
data =
model.generate
({x, y, z}, 20000)
40
41
# Project pdf and data on x
42
# -------------------------------------------------
43
44
# Make plain projection of data and pdf on x observable
45
frame =
x.frame
(Title=
"Projection of 3D data and pdf on X"
, Bins=40)
46
data.plotOn
(frame)
47
model.plotOn
(frame)
48
49
# Project pdf and data on x in signal range
50
# ----------------------------------------------------------------------------------
51
52
# Define signal region in y and z observables
53
y.setRange
(
"sigRegion"
, -1, 1)
54
z.setRange
(
"sigRegion"
, -1, 1)
55
56
# Make plot frame
57
frame2 =
x.frame
(Title=
"Same projection on X in signal range of (Y,Z)"
, Bins=40)
58
59
# Plot subset of data in which all observables are inside "sigRegion"
60
# For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
61
# an implicit definition is used that is identical to the full range (i.e.
62
# [-5,5] for x)
63
data.plotOn
(frame2, CutRange=
"sigRegion"
)
64
65
# Project model on x, projected observables (y,z) only in "sigRegion"
66
model.plotOn
(frame2, ProjectionRange=
"sigRegion"
)
67
68
c =
ROOT.TCanvas
(
"rf311_rangeplot"
,
"rf310_rangeplot"
, 800, 400)
69
c.Divide
(2)
70
c.cd
(1)
71
ROOT.gPad.SetLeftMargin
(0.15)
72
frame.GetYaxis
().SetTitleOffset(1.4)
73
frame.Draw
()
74
c.cd
(2)
75
ROOT.gPad.SetLeftMargin
(0.15)
76
frame2.GetYaxis
().SetTitleOffset(1.4)
77
frame2.Draw
()
78
79
c.SaveAs
(
"rf311_rangeplot.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
rf311_rangeplot.py
ROOT master - Reference Guide Generated on Fri Jan 24 2025 14:35:49 (GVA Time) using Doxygen 1.10.0