ROOT
Version v6.34
master
v6.32
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf301_composition.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Multidimensional models: multi-dimensional pdfs through composition, e.g. substituting
5
## a pdf parameter with a function that depends on other observables
6
##
7
## `pdf = gauss(x,f(y),s)` with `f(y) = a0 + a1*y`
8
##
9
## \macro_image
10
## \macro_code
11
## \macro_output
12
##
13
## \date February 2018
14
## \authors Clemens Lange, Wouter Verkerke (C++ version)
15
16
import
ROOT
17
18
# Setup composed model gauss(x, m(y), s)
19
# -----------------------------------------------------------------------
20
21
# Create observables
22
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -5, 5)
23
y =
ROOT.RooRealVar
(
"y"
,
"y"
, -5, 5)
24
25
# Create function f(y) = a0 + a1*y
26
a0 =
ROOT.RooRealVar
(
"a0"
,
"a0"
, -0.5, -5, 5)
27
a1 =
ROOT.RooRealVar
(
"a1"
,
"a1"
, -0.5, -1, 1)
28
fy =
ROOT.RooPolyVar
(
"fy"
,
"fy"
, y, [a0, a1])
29
30
# Creat gauss(x,f(y),s)
31
sigma =
ROOT.RooRealVar
(
"sigma"
,
"width of gaussian"
, 0.5)
32
model =
ROOT.RooGaussian
(
"model"
,
"Gaussian with shifting mean"
, x, fy, sigma)
33
34
# Sample data, plot data and pdf on x and y
35
# ---------------------------------------------------------------------------------
36
37
# Generate 10000 events in x and y from model
38
data =
model.generate
({x, y}, 10000)
39
40
# Plot x distribution of data and projection of model x = Int(dy)
41
# model(x,y)
42
xframe =
x.frame
()
43
data.plotOn
(xframe)
44
model.plotOn
(xframe)
45
46
# Plot x distribution of data and projection of model y = Int(dx)
47
# model(x,y)
48
yframe =
y.frame
()
49
data.plotOn
(yframe)
50
model.plotOn
(yframe)
51
52
# Make two-dimensional plot in x vs y
53
hh_model =
model.createHistogram
(
"hh_model"
, x, Binning=50, YVar=dict(var=y, Binning=50))
54
hh_model.SetLineColor
(
ROOT.kBlue
)
55
56
# Make canvas and draw ROOT.RooPlots
57
c =
ROOT.TCanvas
(
"rf301_composition"
,
"rf301_composition"
, 1200, 400)
58
c.Divide
(3)
59
c.cd
(1)
60
ROOT.gPad.SetLeftMargin
(0.15)
61
xframe.GetYaxis
().SetTitleOffset(1.4)
62
xframe.Draw
()
63
c.cd
(2)
64
ROOT.gPad.SetLeftMargin
(0.15)
65
yframe.GetYaxis
().SetTitleOffset(1.4)
66
yframe.Draw
()
67
c.cd
(3)
68
ROOT.gPad.SetLeftMargin
(0.20)
69
hh_model.GetZaxis
().SetTitleOffset(2.5)
70
hh_model.Draw
(
"surf"
)
71
72
c.SaveAs
(
"rf301_composition.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
rf301_composition.py
ROOT tags/6-34-04 - Reference Guide Generated on Fri Mar 21 2025 04:40:18 (GVA Time) using Doxygen 1.10.0