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
rf304_uncorrprod.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Multidimensional models: simple uncorrelated multi-dimensional pdfs
5
##
6
## `pdf = gauss(x,mx,sx) * gauss(y,my,sy)`
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
18
# Create component pdfs in x and y
19
# ----------------------------------------------------------------
20
21
# Create two pdfs gaussx(x,meanx,sigmax) gaussy(y,meany,sigmay) and its
22
# variables
23
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -5, 5)
24
y =
ROOT.RooRealVar
(
"y"
,
"y"
, -5, 5)
25
26
meanx =
ROOT.RooRealVar
(
"mean1"
,
"mean of gaussian x"
, 2)
27
meany =
ROOT.RooRealVar
(
"mean2"
,
"mean of gaussian y"
, -2)
28
sigmax =
ROOT.RooRealVar
(
"sigmax"
,
"width of gaussian x"
, 1)
29
sigmay =
ROOT.RooRealVar
(
"sigmay"
,
"width of gaussian y"
, 5)
30
31
gaussx =
ROOT.RooGaussian
(
"gaussx"
,
"gaussian PDF"
, x, meanx, sigmax)
32
gaussy =
ROOT.RooGaussian
(
"gaussy"
,
"gaussian PDF"
, y, meany, sigmay)
33
34
# Construct uncorrelated product pdf
35
# -------------------------------------------------------------------
36
37
# Multiply gaussx and gaussy into a two-dimensional pdf gaussxy
38
gaussxy =
ROOT.RooProdPdf
(
"gaussxy"
,
"gaussx*gaussy"
, [gaussx, gaussy])
39
40
# Sample pdf, plot projection on x and y
41
# ---------------------------------------------------------------------------
42
43
# Generate 10000 events in x and y from gaussxy
44
data =
gaussxy.generate
({x, y}, 10000)
45
46
# Plot x distribution of data and projection of gaussxy x = Int(dy)
47
# gaussxy(x,y)
48
xframe =
x.frame
(Title=
"X projection of gauss(x)*gauss(y)"
)
49
data.plotOn
(xframe)
50
gaussxy.plotOn
(xframe)
51
52
# Plot x distribution of data and projection of gaussxy y = Int(dx)
53
# gaussxy(x,y)
54
yframe =
y.frame
(Title=
"Y projection of gauss(x)*gauss(y)"
)
55
data.plotOn
(yframe)
56
gaussxy.plotOn
(yframe)
57
58
# Make canvas and draw ROOT.RooPlots
59
c =
ROOT.TCanvas
(
"rf304_uncorrprod"
,
"rf304_uncorrprod"
, 800, 400)
60
c.Divide
(2)
61
c.cd
(1)
62
ROOT.gPad.SetLeftMargin
(0.15)
63
xframe.GetYaxis
().SetTitleOffset(1.4)
64
xframe.Draw
()
65
c.cd
(2)
66
ROOT.gPad.SetLeftMargin
(0.15)
67
yframe.GetYaxis
().SetTitleOffset(1.4)
68
yframe.Draw
()
69
70
c.SaveAs
(
"rf304_uncorrprod.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
rf304_uncorrprod.py
ROOT tags/6-34-04 - Reference Guide Generated on Fri Mar 21 2025 04:40:18 (GVA Time) using Doxygen 1.10.0