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
5 ## composite (multidimensional) pdfs
6 ##
7 ## \macro_code
8 ##
9 ## \date February 2018
10 ## \authors Clemens Lange, Wouter Verkerke (C++ version)
11 
12 import ROOT
13 
14 # Create observables, parameters
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", -1.5, -5, 5)
23 a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
24 sigma = 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))
30 fy_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)
34 model_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
41 fy_2 = ROOT.RooPolyVar("fy_2", "fy_2", y, ROOT.RooArgList(a0, a1))
42 
43 # Create gauss(x,f(y),s)
44 model_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
51 fy_3 = ROOT.RooAddition("fy_3", "a0+y", ROOT.RooArgList(a0, y))
52 
53 # Create gauss(x,f(y),s)
54 model_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
61 fy_4 = ROOT.RooProduct("fy_4", "a1*y", ROOT.RooArgList(a1, y))
62 
63 # Create gauss(x,f(y),s)
64 model_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
71 hh_model_1 = model_1.createHistogram("hh_model_1", x, ROOT.RooFit.Binning(
72  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
73 hh_model_2 = model_2.createHistogram("hh_model_2", x, ROOT.RooFit.Binning(
74  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
75 hh_model_3 = model_3.createHistogram("hh_model_3", x, ROOT.RooFit.Binning(
76  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
77 hh_model_4 = model_4.createHistogram("hh_model_4", x, ROOT.RooFit.Binning(
78  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
79 hh_model_1.SetLineColor(ROOT.kBlue)
80 hh_model_2.SetLineColor(ROOT.kBlue)
81 hh_model_3.SetLineColor(ROOT.kBlue)
82 hh_model_4.SetLineColor(ROOT.kBlue)
83 
84 # Make canvas and draw ROOT.RooPlots
85 c = ROOT.TCanvas("rf302_utilfuncs", "rf302_utilfuncs", 800, 800)
86 c.Divide(2, 2)
87 c.cd(1)
88 ROOT.gPad.SetLeftMargin(0.20)
89 hh_model_1.GetZaxis().SetTitleOffset(2.5)
90 hh_model_1.Draw("surf")
91 c.cd(2)
92 ROOT.gPad.SetLeftMargin(0.20)
93 hh_model_2.GetZaxis().SetTitleOffset(2.5)
94 hh_model_2.Draw("surf")
95 c.cd(3)
96 ROOT.gPad.SetLeftMargin(0.20)
97 hh_model_3.GetZaxis().SetTitleOffset(2.5)
98 hh_model_3.Draw("surf")
99 c.cd(4)
100 ROOT.gPad.SetLeftMargin(0.20)
101 hh_model_4.GetZaxis().SetTitleOffset(2.5)
102 hh_model_4.Draw("surf")
103 
104 c.SaveAs("rf302_utilfuncs.png")