Logo ROOT   6.16/01
Reference Guide
rf309_ndimplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #308
6##
7## Making 2/3 dimensional plots of p.d.f.s and datasets
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15import ROOT
16
17# Create 2D model and dataset
18# -----------------------------------------------------
19
20# Create observables
21x = ROOT.RooRealVar("x", "x", -5, 5)
22y = ROOT.RooRealVar("y", "y", -5, 5)
23
24# Create parameters
25a0 = ROOT.RooRealVar("a0", "a0", -3.5, -5, 5)
26a1 = ROOT.RooRealVar("a1", "a1", -1.5, -1, 1)
27sigma = ROOT.RooRealVar("sigma", "width of gaussian", 1.5)
28
29# Create interpreted function f(y) = a0 - a1*sqrt(10*abs(y))
30fy = ROOT.RooFormulaVar("fy", "a0-a1*sqrt(10*abs(y))",
31 ROOT.RooArgList(y, a0, a1))
32
33# Create gauss(x,f(y),s)
34model = ROOT.RooGaussian(
35 "model", "Gaussian with shifting mean", x, fy, sigma)
36
37# Sample dataset from gauss(x,y)
38data = model.generate(ROOT.RooArgSet(x, y), 10000)
39
40# Make 2D plots of data and model
41# -------------------------------------------------------------
42
43# Create and fill ROOT 2D histogram (20x20 bins) with contents of dataset
44# hh_data = data.createHistogram("hh_data",x, ROOT.RooFit.Binning(20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
45# hh_data = data.createHistogram("x,y", 20, 20) # does not work, see
46# https://root.cern.ch/phpBB3/viewtopic.php?t=16648
47hh_data = ROOT.RooAbsData.createHistogram(data, "x,y", x, ROOT.RooFit.Binning(
48 20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
49
50# Create and fill ROOT 2D histogram (50x50 bins) with sampling of pdf
51# hh_pdf = model.createHistogram("hh_model",x, ROOT.RooFit.Binning(50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
52hh_pdf = model.createHistogram("x,y", 50, 50)
53hh_pdf.SetLineColor(ROOT.kBlue)
54
55# Create 3D model and dataset
56# -----------------------------------------------------
57
58# Create observables
59z = ROOT.RooRealVar("z", "z", -5, 5)
60
61gz = ROOT.RooGaussian(
62 "gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))
63model3 = ROOT.RooProdPdf("model3", "model3", ROOT.RooArgList(model, gz))
64
65data3 = model3.generate(ROOT.RooArgSet(x, y, z), 10000)
66
67# Make 3D plots of data and model
68# -------------------------------------------------------------
69
70# Create and fill ROOT 2D histogram (8x8x8 bins) with contents of dataset
71# hh_data3 = data3.createHistogram("hh_data3", x, ROOT.RooFit.Binning(8), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(8)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(8)))
72hh_data3 = ROOT.RooAbsData.createHistogram(data3, "hh_data3", x, ROOT.RooFit.Binning(
73 8), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(8)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(8)))
74
75# Create and fill ROOT 2D histogram (20x20x20 bins) with sampling of pdf
76hh_pdf3 = model3.createHistogram("hh_model3", x, ROOT.RooFit.Binning(
77 20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(20)))
78hh_pdf3.SetFillColor(ROOT.kBlue)
79
80c1 = ROOT.TCanvas("rf309_2dimplot", "rf309_2dimplot", 800, 800)
81c1.Divide(2, 2)
82c1.cd(1)
83ROOT.gPad.SetLeftMargin(0.15)
84hh_data.GetZaxis().SetTitleOffset(1.4)
85hh_data.Draw("lego")
86c1.cd(2)
87ROOT.gPad.SetLeftMargin(0.20)
88hh_pdf.GetZaxis().SetTitleOffset(2.5)
89hh_pdf.Draw("surf")
90c1.cd(3)
91ROOT.gPad.SetLeftMargin(0.15)
92hh_data.GetZaxis().SetTitleOffset(1.4)
93hh_data.Draw("box")
94c1.cd(4)
95ROOT.gPad.SetLeftMargin(0.15)
96hh_pdf.GetZaxis().SetTitleOffset(2.5)
97hh_pdf.Draw("cont3")
98c1.SaveAs("rf309_2dimplot.png")
99
100c2 = ROOT.TCanvas("rf309_3dimplot", "rf309_3dimplot", 800, 400)
101c2.Divide(2)
102c2.cd(1)
103ROOT.gPad.SetLeftMargin(0.15)
104hh_data3.GetZaxis().SetTitleOffset(1.4)
105hh_data3.Draw("lego")
106c2.cd(2)
107ROOT.gPad.SetLeftMargin(0.15)
108hh_pdf3.GetZaxis().SetTitleOffset(1.4)
109hh_pdf3.Draw("iso")
110c2.SaveAs("rf309_3dimplot.png")