Logo ROOT  
Reference Guide
rf103_interprfuncs.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Basic functionality: interpreted functions and p.d.f.s
5##
6## \macro_code
7##
8## \date February 2018
9## \author Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Generic interpreted p.d.f.
14# ------------------------------
15
16# Declare observable x
17x = ROOT.RooRealVar("x", "x", -20, 20)
18
19# Construct generic pdf from interpreted expression
20# ------------------------------------------------------
21
22# ROOT.To construct a proper p.d.f, the formula expression is explicitly normalized internally by dividing
23# it by a numeric integral of the expresssion over x in the range [-20,20]
24#
25alpha = ROOT.RooRealVar("alpha", "alpha", 5, 0.1, 10)
26genpdf = ROOT.RooGenericPdf(
27 "genpdf",
28 "genpdf",
29 "(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))",
30 ROOT.RooArgList(
31 x,
32 alpha))
33
34# Sample, fit and plot generic pdf
35# ---------------------------------------------------------------
36
37# Generate a toy dataset from the interpreted p.d.f
38data = genpdf.generate(ROOT.RooArgSet(x), 10000)
39
40# Fit the interpreted p.d.f to the generated data
41genpdf.fitTo(data)
42
43# Make a plot of the data and the p.d.f overlaid
44xframe = x.frame(ROOT.RooFit.Title("Interpreted expression pdf"))
45data.plotOn(xframe)
46genpdf.plotOn(xframe)
47
48# Standard p.d.f. adjust with interpreted helper function
49# ------------------------------------------------------------------------------------------------------------
50# Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian #
51#
52# Construct standard pdf with formula replacing parameter
53# ------------------------------------------------------------------------------------------------------------
54
55# Construct parameter mean2 and sigma
56mean2 = ROOT.RooRealVar("mean2", "mean^2", 10, 0, 200)
57sigma = ROOT.RooRealVar("sigma", "sigma", 3, 0.1, 10)
58
59# Construct interpreted function mean = sqrt(mean^2)
60mean = ROOT.RooFormulaVar(
61 "mean", "mean", "sqrt(mean2)", ROOT.RooArgList(mean2))
62
63# Construct a gaussian g2(x,sqrt(mean2),sigma)
64g2 = ROOT.RooGaussian("g2", "h2", x, mean, sigma)
65
66# Generate toy data
67# ---------------------------------
68
69# Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian
70# dataset with mean 10 and width 3
71g1 = ROOT.RooGaussian("g1", "g1", x, ROOT.RooFit.RooConst(
72 10), ROOT.RooFit.RooConst(3))
73data2 = g1.generate(ROOT.RooArgSet(x), 1000)
74
75# Fit and plot tailored standard pdf
76# -------------------------------------------------------------------
77
78# Fit g2 to data from g1
79r = g2.fitTo(data2, ROOT.RooFit.Save()) # ROOT.RooFitResult
80r.Print()
81
82# Plot data on frame and overlay projection of g2
83xframe2 = x.frame(ROOT.RooFit.Title("Tailored Gaussian pdf"))
84data2.plotOn(xframe2)
85g2.plotOn(xframe2)
86
87# Draw all frames on a canvas
88c = ROOT.TCanvas("rf103_interprfuncs", "rf103_interprfuncs", 800, 400)
89c.Divide(2)
90c.cd(1)
91ROOT.gPad.SetLeftMargin(0.15)
92xframe.GetYaxis().SetTitleOffset(1.4)
93xframe.Draw()
94c.cd(2)
95ROOT.gPad.SetLeftMargin(0.15)
96xframe2.GetYaxis().SetTitleOffset(1.4)
97xframe2.Draw()
98
99c.SaveAs("rf103_interprfuncs.png")