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