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