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("fy_1", "a0-a1*sqrt(10*abs(y))", ROOT.RooArgList(y, a0, a1))
31 
32 # Create gauss(x,f(y),s)
33 model_1 = ROOT.RooGaussian("model_1", "Gaussian with shifting mean", x, fy_1, sigma)
34 
35 # Using RooPolyVar to tailor pdf
36 # -----------------------------------------------------------------------
37 
38 # Create polynomial function f(y) = a0 + a1*y
39 fy_2 = ROOT.RooPolyVar("fy_2", "fy_2", y, ROOT.RooArgList(a0, a1))
40 
41 # Create gauss(x,f(y),s)
42 model_2 = ROOT.RooGaussian("model_2", "Gaussian with shifting mean", x, fy_2, sigma)
43 
44 # Using RooAddition to tailor pdf
45 # -----------------------------------------------------------------------
46 
47 # Create sum function f(y) = a0 + y
48 fy_3 = ROOT.RooAddition("fy_3", "a0+y", ROOT.RooArgList(a0, y))
49 
50 # Create gauss(x,f(y),s)
51 model_3 = ROOT.RooGaussian("model_3", "Gaussian with shifting mean", x, fy_3, sigma)
52 
53 # Using RooProduct to tailor pdf
54 # -----------------------------------------------------------------------
55 
56 # Create product function f(y) = a1*y
57 fy_4 = ROOT.RooProduct("fy_4", "a1*y", ROOT.RooArgList(a1, y))
58 
59 # Create gauss(x,f(y),s)
60 model_4 = ROOT.RooGaussian("model_4", "Gaussian with shifting mean", x, fy_4, sigma)
61 
62 # Plot all pdfs
63 # ----------------------------
64 
65 # Make two-dimensional plots in x vs y
66 hh_model_1 = model_1.createHistogram("hh_model_1", x, Binning=50, YVar=(y, ROOT.RooFit.Binning(50)))
67 hh_model_2 = model_2.createHistogram("hh_model_2", x, Binning=50, YVar=(y, ROOT.RooFit.Binning(50)))
68 hh_model_3 = model_3.createHistogram("hh_model_3", x, Binning=50, YVar=(y, ROOT.RooFit.Binning(50)))
69 hh_model_4 = model_4.createHistogram("hh_model_4", x, Binning=50, YVar=(y, ROOT.RooFit.Binning(50)))
70 hh_model_1.SetLineColor(ROOT.kBlue)
71 hh_model_2.SetLineColor(ROOT.kBlue)
72 hh_model_3.SetLineColor(ROOT.kBlue)
73 hh_model_4.SetLineColor(ROOT.kBlue)
74 
75 # Make canvas and draw ROOT.RooPlots
76 c = ROOT.TCanvas("rf302_utilfuncs", "rf302_utilfuncs", 800, 800)
77 c.Divide(2, 2)
78 c.cd(1)
79 ROOT.gPad.SetLeftMargin(0.20)
80 hh_model_1.GetZaxis().SetTitleOffset(2.5)
81 hh_model_1.Draw("surf")
82 c.cd(2)
83 ROOT.gPad.SetLeftMargin(0.20)
84 hh_model_2.GetZaxis().SetTitleOffset(2.5)
85 hh_model_2.Draw("surf")
86 c.cd(3)
87 ROOT.gPad.SetLeftMargin(0.20)
88 hh_model_3.GetZaxis().SetTitleOffset(2.5)
89 hh_model_3.Draw("surf")
90 c.cd(4)
91 ROOT.gPad.SetLeftMargin(0.20)
92 hh_model_4.GetZaxis().SetTitleOffset(2.5)
93 hh_model_4.Draw("surf")
94 
95 c.SaveAs("rf302_utilfuncs.png")