Logo ROOT  
Reference Guide
 
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_code
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Create component pdfs in x and y
17# ----------------------------------------------------------------
18
19# Create two pdfs 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 pdf gaussxy
36gaussxy = ROOT.RooProdPdf("gaussxy", "gaussx*gaussy", [gaussx, gaussy])
37
38# Sample pdf, plot projection on x and y
39# ---------------------------------------------------------------------------
40
41# Generate 10000 events in x and y from gaussxy
42data = gaussxy.generate({x, y}, 10000)
43
44# Plot x distribution of data and projection of gaussxy x = Int(dy)
45# gaussxy(x,y)
46xframe = x.frame(Title="X projection of gauss(x)*gauss(y)")
47data.plotOn(xframe)
48gaussxy.plotOn(xframe)
49
50# Plot x distribution of data and projection of gaussxy y = Int(dx)
51# gaussxy(x,y)
52yframe = y.frame(Title="Y projection of gauss(x)*gauss(y)")
53data.plotOn(yframe)
54gaussxy.plotOn(yframe)
55
56# Make canvas and draw ROOT.RooPlots
57c = ROOT.TCanvas("rf304_uncorrprod", "rf304_uncorrprod", 800, 400)
58c.Divide(2)
59c.cd(1)
60ROOT.gPad.SetLeftMargin(0.15)
61xframe.GetYaxis().SetTitleOffset(1.4)
62xframe.Draw()
63c.cd(2)
64ROOT.gPad.SetLeftMargin(0.15)
65yframe.GetYaxis().SetTitleOffset(1.4)
66yframe.Draw()
67
68c.SaveAs("rf304_uncorrprod.png")