Logo ROOT  
Reference Guide
rf302_utilfuncs.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: utility functions classes available for use in tailoring of composite (multidimensional) pdfs
5##
6## \macro_code
7##
8## \date February 2018
9## \author Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Create observables, parameters
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", -1.5, -5, 5)
22a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
23sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
24
25# Using RooFormulaVar to tailor pdf
26# -----------------------------------------------------------------------
27
28# Create interpreted function f(y) = a0 - a1*sqrt(10*abs(y))
29fy_1 = ROOT.RooFormulaVar(
30 "fy_1", "a0-a1*sqrt(10*abs(y))", ROOT.RooArgList(y, a0, a1))
31
32# Create gauss(x,f(y),s)
33model_1 = ROOT.RooGaussian(
34 "model_1", "Gaussian with shifting mean", x, fy_1, sigma)
35
36# Using RooPolyVar to tailor pdf
37# -----------------------------------------------------------------------
38
39# Create polynomial function f(y) = a0 + a1*y
40fy_2 = ROOT.RooPolyVar("fy_2", "fy_2", y, ROOT.RooArgList(a0, a1))
41
42# Create gauss(x,f(y),s)
43model_2 = ROOT.RooGaussian(
44 "model_2", "Gaussian with shifting mean", x, fy_2, sigma)
45
46# Using RooAddition to tailor pdf
47# -----------------------------------------------------------------------
48
49# Create sum function f(y) = a0 + y
50fy_3 = ROOT.RooAddition("fy_3", "a0+y", ROOT.RooArgList(a0, y))
51
52# Create gauss(x,f(y),s)
53model_3 = ROOT.RooGaussian(
54 "model_3", "Gaussian with shifting mean", x, fy_3, sigma)
55
56# Using RooProduct to tailor pdf
57# -----------------------------------------------------------------------
58
59# Create product function f(y) = a1*y
60fy_4 = ROOT.RooProduct("fy_4", "a1*y", ROOT.RooArgList(a1, y))
61
62# Create gauss(x,f(y),s)
63model_4 = ROOT.RooGaussian(
64 "model_4", "Gaussian with shifting mean", x, fy_4, sigma)
65
66# Plot all pdfs
67# ----------------------------
68
69# Make two-dimensional plots in x vs y
70hh_model_1 = model_1.createHistogram("hh_model_1", x, ROOT.RooFit.Binning(
71 50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
72hh_model_2 = model_2.createHistogram("hh_model_2", x, ROOT.RooFit.Binning(
73 50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
74hh_model_3 = model_3.createHistogram("hh_model_3", x, ROOT.RooFit.Binning(
75 50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
76hh_model_4 = model_4.createHistogram("hh_model_4", x, ROOT.RooFit.Binning(
77 50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
78hh_model_1.SetLineColor(ROOT.kBlue)
79hh_model_2.SetLineColor(ROOT.kBlue)
80hh_model_3.SetLineColor(ROOT.kBlue)
81hh_model_4.SetLineColor(ROOT.kBlue)
82
83# Make canvas and draw ROOT.RooPlots
84c = ROOT.TCanvas("rf302_utilfuncs", "rf302_utilfuncs", 800, 800)
85c.Divide(2, 2)
86c.cd(1)
87ROOT.gPad.SetLeftMargin(0.20)
88hh_model_1.GetZaxis().SetTitleOffset(2.5)
89hh_model_1.Draw("surf")
90c.cd(2)
91ROOT.gPad.SetLeftMargin(0.20)
92hh_model_2.GetZaxis().SetTitleOffset(2.5)
93hh_model_2.Draw("surf")
94c.cd(3)
95ROOT.gPad.SetLeftMargin(0.20)
96hh_model_3.GetZaxis().SetTitleOffset(2.5)
97hh_model_3.Draw("surf")
98c.cd(4)
99ROOT.gPad.SetLeftMargin(0.20)
100hh_model_4.GetZaxis().SetTitleOffset(2.5)
101hh_model_4.Draw("surf")
102
103c.SaveAs("rf302_utilfuncs.png")