ROOT
Version master
v6.34
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
rf305_condcorrprod.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
3
## \notebook
4
## Multidimensional models: multi-dimensional pdfs with conditional pdfs in product
5
##
6
## `pdf = gauss(x,f(y),sx | y ) * gauss(y,ms,sx)` with `f(y) = a0 + a1*y`
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C++ version)
14
15
import
ROOT
16
17
# Create conditional pdf gx(x|y)
18
# -----------------------------------------------------------
19
20
# Create observables
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -5, 5)
22
y =
ROOT.RooRealVar
(
"y"
,
"y"
, -5, 5)
23
24
# Create function f(y) = a0 + a1*y
25
a0 =
ROOT.RooRealVar
(
"a0"
,
"a0"
, -0.5, -5, 5)
26
a1 =
ROOT.RooRealVar
(
"a1"
,
"a1"
, -0.5, -1, 1)
27
fy =
ROOT.RooPolyVar
(
"fy"
,
"fy"
, y, [a0, a1])
28
29
# Create gaussx(x,f(y),sx)
30
sigmax =
ROOT.RooRealVar
(
"sigma"
,
"width of gaussian"
, 0.5)
31
gaussx =
ROOT.RooGaussian
(
"gaussx"
,
"Gaussian in x with shifting mean in y"
, x, fy, sigmax)
32
33
# Create pdf gy(y)
34
# -----------------------------------------------------------
35
36
# Create gaussy(y,0,5)
37
gaussy =
ROOT.RooGaussian
(
"gaussy"
,
"Gaussian in y"
, y, 0.0, 3.0)
38
39
# Create product gx(x|y)*gy(y)
40
# -------------------------------------------------------
41
42
# Create gaussx(x,sx|y) * gaussy(y)
43
model =
ROOT.RooProdPdf
(
"model"
,
"gaussx(x|y)*gaussy(y)"
, {gaussy}, Conditional=({gaussx}, {x}))
44
45
# Sample, fit and plot product pdf
46
# ---------------------------------------------------------------
47
48
# Generate 1000 events in x and y from model
49
data =
model.generate
({x, y}, 10000)
50
51
# Plot x distribution of data and projection of model x = Int(dy)
52
# model(x,y)
53
xframe =
x.frame
()
54
data.plotOn
(xframe)
55
model.plotOn
(xframe)
56
57
# Plot x distribution of data and projection of model y = Int(dx)
58
# model(x,y)
59
yframe =
y.frame
()
60
data.plotOn
(yframe)
61
model.plotOn
(yframe)
62
63
# Make two-dimensional plot in x vs y
64
hh_model =
model.createHistogram
(
"hh_model"
, x,
ROOT.RooFit.Binning
(50),
ROOT.RooFit.YVar
(y,
ROOT.RooFit.Binning
(50)))
65
hh_model.SetLineColor
(
"kBlue"
)
66
67
# Make canvas and draw ROOT.RooPlots
68
c =
ROOT.TCanvas
(
"rf305_condcorrprod"
,
"rf05_condcorrprod"
, 1200, 400)
69
c.Divide
(3)
70
c.cd
(1)
71
ROOT.gPad.SetLeftMargin
(0.15)
72
xframe.GetYaxis
().SetTitleOffset(1.6)
73
xframe.Draw
()
74
c.cd
(2)
75
ROOT.gPad.SetLeftMargin
(0.15)
76
yframe.GetYaxis
().SetTitleOffset(1.6)
77
yframe.Draw
()
78
c.cd
(3)
79
ROOT.gPad.SetLeftMargin
(0.20)
80
hh_model.GetZaxis
().SetTitleOffset(2.5)
81
hh_model.Draw
(
"surf"
)
82
83
c.SaveAs
(
"rf305_condcorrprod.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
rf305_condcorrprod.py
ROOT master - Reference Guide Generated on Tue Apr 8 2025 06:21:58 (GVA Time) using Doxygen 1.10.0