rf304_uncorrprod

Multidimensional models: simple uncorrelated multi-dimensional pdfs

pdf = gauss(x,mx,sx) * gauss(y,my,sy)

Author: Clemens Lange, Wouter Verkerke (C++ version)
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, November 30, 2022 at 11:22 AM.

In [1]:
import ROOT
Welcome to JupyROOT 6.27/01

Create component pdfs in x and y

Create two pdfs gaussx(x,meanx,sigmax) gaussy(y,meany,sigmay) and its variables

In [2]:
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -5, 5)

meanx = ROOT.RooRealVar("mean1", "mean of gaussian x", 2)
meany = ROOT.RooRealVar("mean2", "mean of gaussian y", -2)
sigmax = ROOT.RooRealVar("sigmax", "width of gaussian x", 1)
sigmay = ROOT.RooRealVar("sigmay", "width of gaussian y", 5)

gaussx = ROOT.RooGaussian("gaussx", "gaussian PDF", x, meanx, sigmax)
gaussy = ROOT.RooGaussian("gaussy", "gaussian PDF", y, meany, sigmay)
[#0] WARNING:InputArguments -- The parameter 'sigmax' with range [-1e+30, 1e+30] of the RooGaussian 'gaussx' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigmay' with range [-1e+30, 1e+30] of the RooGaussian 'gaussy' exceeds the safe range of (0, inf). Advise to limit its range.

Construct uncorrelated product pdf

Multiply gaussx and gaussy into a two-dimensional pdf gaussxy

In [3]:
gaussxy = ROOT.RooProdPdf("gaussxy", "gaussx*gaussy", [gaussx, gaussy])

Sample pdf, plot projection on x and y

Generate 10000 events in x and y from gaussxy

In [4]:
data = gaussxy.generate({x, y}, 10000)

Plot x distribution of data and projection of gaussxy x = Int(dy) gaussxy(x,y)

In [5]:
xframe = x.frame(Title="X projection of gauss(x)*gauss(y)")
data.plotOn(xframe)
gaussxy.plotOn(xframe)
Out[5]:
<cppyy.gbl.RooPlot object at 0x94d32c0>
[#1] INFO:Plotting -- RooAbsReal::plotOn(gaussxy) plot on x integrates over variables (y)

Plot x distribution of data and projection of gaussxy y = Int(dx) gaussxy(x,y)

In [6]:
yframe = y.frame(Title="Y projection of gauss(x)*gauss(y)")
data.plotOn(yframe)
gaussxy.plotOn(yframe)
Out[6]:
<cppyy.gbl.RooPlot object at 0x96a9d90>
[#1] INFO:Plotting -- RooAbsReal::plotOn(gaussxy) plot on y integrates over variables (x)

Make canvas and draw ROOT.RooPlots

In [7]:
c = ROOT.TCanvas("rf304_uncorrprod", "rf304_uncorrprod", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
yframe.GetYaxis().SetTitleOffset(1.4)
yframe.Draw()

c.SaveAs("rf304_uncorrprod.png")
Info in <TCanvas::Print>: png file rf304_uncorrprod.png has been created

Draw all canvases

In [8]:
from ROOT import gROOT 
gROOT.GetListOfCanvases().Draw()