Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf309_ndimplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: making 2/3 dimensional plots of pdfs and datasets
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Create 2D model and dataset
14# -----------------------------------------------------
15
16# Create observables
17x = ROOT.RooRealVar("x", "x", -5, 5)
18y = ROOT.RooRealVar("y", "y", -5, 5)
19
20# Create parameters
21a0 = ROOT.RooRealVar("a0", "a0", -3.5, -5, 5)
22a1 = ROOT.RooRealVar("a1", "a1", -1.5, -1, 1)
23sigma = ROOT.RooRealVar("sigma", "width of gaussian", 1.5)
24
25# Create interpreted function f(y) = a0 - a1*sqrt(10*abs(y))
26fy = ROOT.RooFormulaVar("fy", "a0-a1*sqrt(10*abs(y))",
27 ROOT.RooArgList(y, a0, a1))
28
29# Create gauss(x,f(y),s)
30model = ROOT.RooGaussian(
31 "model", "Gaussian with shifting mean", x, fy, sigma)
32
33# Sample dataset from gauss(x,y)
34data = model.generate(ROOT.RooArgSet(x, y), 10000)
35
36# Make 2D plots of data and model
37# -------------------------------------------------------------
38
39# Create and fill ROOT 2D histogram (20x20 bins) with contents of dataset
40# hh_data = data.createHistogram("hh_data",x, ROOT.RooFit.Binning(20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
41# hh_data = data.createHistogram("x,y", 20, 20) # does not work, see
42# https://root.cern.ch/phpBB3/viewtopic.php?t=16648
43hh_data = ROOT.RooAbsData.createHistogram(data, "x,y", x, ROOT.RooFit.Binning(
44 20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
45
46# Create and fill ROOT 2D histogram (50x50 bins) with sampling of pdf
47# hh_pdf = model.createHistogram("hh_model",x, ROOT.RooFit.Binning(50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
48hh_pdf = model.createHistogram("x,y", 50, 50)
49hh_pdf.SetLineColor(ROOT.kBlue)
50
51# Create 3D model and dataset
52# -----------------------------------------------------
53
54# Create observables
55z = ROOT.RooRealVar("z", "z", -5, 5)
56
57gz = ROOT.RooGaussian(
58 "gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))
59model3 = ROOT.RooProdPdf("model3", "model3", ROOT.RooArgList(model, gz))
60
61data3 = model3.generate(ROOT.RooArgSet(x, y, z), 10000)
62
63# Make 3D plots of data and model
64# -------------------------------------------------------------
65
66# Create and fill ROOT 2D histogram (8x8x8 bins) with contents of dataset
67# 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)))
68hh_data3 = ROOT.RooAbsData.createHistogram(
69 data3, "hh_data3", x, ROOT.RooFit.Binning(8), ROOT.RooFit.YVar(
70 y, ROOT.RooFit.Binning(8)), ROOT.RooFit.ZVar(
71 z, ROOT.RooFit.Binning(8)))
72
73# Create and fill ROOT 2D histogram (20x20x20 bins) with sampling of pdf
74hh_pdf3 = model3.createHistogram(
75 "hh_model3", x, ROOT.RooFit.Binning(20), ROOT.RooFit.YVar(
76 y, ROOT.RooFit.Binning(20)), ROOT.RooFit.ZVar(
77 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")