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))", [y, a0, a1])
27
28# Create gauss(x,f(y),s)
29model = ROOT.RooGaussian("model", "Gaussian with shifting mean", x, fy, sigma)
30
31# Sample dataset from gauss(x,y)
32data = model.generate({x, y}, 10000)
33
34# Make 2D plots of data and model
35# -------------------------------------------------------------
36
37# Create and fill ROOT 2D histogram (20x20 bins) with contents of dataset
38# hh_data = data.createHistogram("hh_data",x, ROOT.RooFit.Binning(20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
39# hh_data = data.createHistogram("x,y", 20, 20) # does not work, see
40# https://root.cern.ch/phpBB3/viewtopic.php?t=16648
41hh_data = ROOT.RooAbsData.createHistogram(data, "x,y", x, Binning=20, YVar=dict(var=y, Binning=20))
42
43# Create and fill ROOT 2D histogram (50x50 bins) with sampling of pdf
44# hh_pdf = model.createHistogram("hh_model",x, ROOT.RooFit.Binning(50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
45hh_pdf = model.createHistogram("x,y", 50, 50)
46hh_pdf.SetLineColor(ROOT.kBlue)
47
48# Create 3D model and dataset
49# -----------------------------------------------------
50
51# Create observables
52z = ROOT.RooRealVar("z", "z", -5, 5)
53
54gz = ROOT.RooGaussian("gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))
55model3 = ROOT.RooProdPdf("model3", "model3", [model, gz])
56
57data3 = model3.generate({x, y, z}, 10000)
58
59# Make 3D plots of data and model
60# -------------------------------------------------------------
61
62# Create and fill ROOT 2D histogram (8x8x8 bins) with contents of dataset
63# 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)))
64hh_data3 = ROOT.RooAbsData.createHistogram(
65 data3,
66 "hh_data3",
67 x,
68 Binning=8,
69 YVar=dict(var=y, Binning=8),
70 ZVar=dict(var=z, Binning=8),
71)
72
73# Create and fill ROOT 2D histogram (20x20x20 bins) with sampling of pdf
74hh_pdf3 = model3.createHistogram(
75 "hh_model3",
76 x,
77 Binning=20,
78 YVar=dict(var=y, Binning=20),
79 ZVar=dict(var=z, Binning=20),
80)
81hh_pdf3.SetFillColor(ROOT.kBlue)
82
83c1 = ROOT.TCanvas("rf309_2dimplot", "rf309_2dimplot", 800, 800)
84c1.Divide(2, 2)
85c1.cd(1)
86ROOT.gPad.SetLeftMargin(0.15)
87hh_data.GetZaxis().SetTitleOffset(1.4)
88hh_data.Draw("lego")
89c1.cd(2)
90ROOT.gPad.SetLeftMargin(0.20)
91hh_pdf.GetZaxis().SetTitleOffset(2.5)
92hh_pdf.Draw("surf")
93c1.cd(3)
94ROOT.gPad.SetLeftMargin(0.15)
95hh_data.GetZaxis().SetTitleOffset(1.4)
96hh_data.Draw("box")
97c1.cd(4)
98ROOT.gPad.SetLeftMargin(0.15)
99hh_pdf.GetZaxis().SetTitleOffset(2.5)
100hh_pdf.Draw("cont3")
101c1.SaveAs("rf309_2dimplot.png")
102
103c2 = ROOT.TCanvas("rf309_3dimplot", "rf309_3dimplot", 800, 400)
104c2.Divide(2)
105c2.cd(1)
106ROOT.gPad.SetLeftMargin(0.15)
107hh_data3.GetZaxis().SetTitleOffset(1.4)
108hh_data3.Draw("lego")
109c2.cd(2)
110ROOT.gPad.SetLeftMargin(0.15)
111hh_pdf3.GetZaxis().SetTitleOffset(1.4)
112hh_pdf3.Draw("iso")
113c2.SaveAs("rf309_3dimplot.png")