Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
13 import 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
21 x = ROOT.RooRealVar("x", "x", -5, 5)
22 y = ROOT.RooRealVar("y", "y", -5, 5)
23 
24 meanx = ROOT.RooRealVar("mean1", "mean of gaussian x", 2)
25 meany = ROOT.RooRealVar("mean2", "mean of gaussian y", -2)
26 sigmax = ROOT.RooRealVar("sigmax", "width of gaussian x", 1)
27 sigmay = ROOT.RooRealVar("sigmay", "width of gaussian y", 5)
28 
29 gaussx = ROOT.RooGaussian("gaussx", "gaussian PDF", x, meanx, sigmax)
30 gaussy = 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
36 gaussxy = 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
43 data = 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)
47 xframe = x.frame(ROOT.RooFit.Title("X projection of gauss(x)*gauss(y)"))
48 data.plotOn(xframe)
49 gaussxy.plotOn(xframe)
50 
51 # Plot x distribution of data and projection of gaussxy y = Int(dx)
52 # gaussxy(x,y)
53 yframe = y.frame(ROOT.RooFit.Title("Y projection of gauss(x)*gauss(y)"))
54 data.plotOn(yframe)
55 gaussxy.plotOn(yframe)
56 
57 # Make canvas and draw ROOT.RooPlots
58 c = ROOT.TCanvas("rf304_uncorrprod", "rf304_uncorrprod", 800, 400)
59 c.Divide(2)
60 c.cd(1)
61 ROOT.gPad.SetLeftMargin(0.15)
62 xframe.GetYaxis().SetTitleOffset(1.4)
63 xframe.Draw()
64 c.cd(2)
65 ROOT.gPad.SetLeftMargin(0.15)
66 yframe.GetYaxis().SetTitleOffset(1.4)
67 yframe.Draw()
68 
69 c.SaveAs("rf304_uncorrprod.png")