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