Logo ROOT   6.18/05
Reference Guide
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 p.d.f.s
5##
6## pdf = gauss(x,mx,sx) * gauss(y,my,sy)
7##
8## \macro_code
9##
10## \date February 2018
11## \author Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Create component pdfs in x and y
17# ----------------------------------------------------------------
18
19# Create two p.d.f.s gaussx(x,meanx,sigmax) gaussy(y,meany,sigmay) and its
20# variables
21x = ROOT.RooRealVar("x", "x", -5, 5)
22y = ROOT.RooRealVar("y", "y", -5, 5)
23
24meanx = ROOT.RooRealVar("mean1", "mean of gaussian x", 2)
25meany = ROOT.RooRealVar("mean2", "mean of gaussian y", -2)
26sigmax = ROOT.RooRealVar("sigmax", "width of gaussian x", 1)
27sigmay = ROOT.RooRealVar("sigmay", "width of gaussian y", 5)
28
29gaussx = ROOT.RooGaussian("gaussx", "gaussian PDF", x, meanx, sigmax)
30gaussy = ROOT.RooGaussian("gaussy", "gaussian PDF", y, meany, sigmay)
31
32# Construct uncorrelated product pdf
33# -------------------------------------------------------------------
34
35# Multiply gaussx and gaussy into a two-dimensional p.d.f. gaussxy
36gaussxy = ROOT.RooProdPdf(
37 "gaussxy", "gaussx*gaussy", ROOT.RooArgList(gaussx, gaussy))
38
39# Sample pdf, plot projection on x and y
40# ---------------------------------------------------------------------------
41
42# Generate 10000 events in x and y from gaussxy
43data = gaussxy.generate(ROOT.RooArgSet(x, y), 10000)
44
45# Plot x distribution of data and projection of gaussxy x = Int(dy)
46# gaussxy(x,y)
47xframe = x.frame(ROOT.RooFit.Title("X projection of gauss(x)*gauss(y)"))
48data.plotOn(xframe)
49gaussxy.plotOn(xframe)
50
51# Plot x distribution of data and projection of gaussxy y = Int(dx)
52# gaussxy(x,y)
53yframe = y.frame(ROOT.RooFit.Title("Y projection of gauss(x)*gauss(y)"))
54data.plotOn(yframe)
55gaussxy.plotOn(yframe)
56
57# Make canvas and draw ROOT.RooPlots
58c = ROOT.TCanvas("rf304_uncorrprod", "rf304_uncorrprod", 800, 400)
59c.Divide(2)
60c.cd(1)
61ROOT.gPad.SetLeftMargin(0.15)
62xframe.GetYaxis().SetTitleOffset(1.4)
63xframe.Draw()
64c.cd(2)
65ROOT.gPad.SetLeftMargin(0.15)
66yframe.GetYaxis().SetTitleOffset(1.4)
67yframe.Draw()
68
69c.SaveAs("rf304_uncorrprod.png")