Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
rf103_interprfuncs.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Basic functionality: interpreted functions and pdfs
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Generic interpreted pdf
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 pdf, 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("genpdf", "genpdf", "(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))", [x, alpha])
27
28# Sample, fit and plot generic pdf
29# ---------------------------------------------------------------
30
31# Generate a toy dataset from the interpreted pdf
32data = genpdf.generate({x}, 10000)
33
34# Fit the interpreted pdf to the generated data
35genpdf.fitTo(data)
36
37# Make a plot of the data and the pdf overlaid
38xframe = x.frame(Title="Interpreted expression pdf")
39data.plotOn(xframe)
40genpdf.plotOn(xframe)
41
42# Standard pdf adjust with interpreted helper function
43# ------------------------------------------------------------------------------------------------------------
44# Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian #
45#
46# Construct standard pdf with formula replacing parameter
47# ------------------------------------------------------------------------------------------------------------
48
49# Construct parameter mean2 and sigma
50mean2 = ROOT.RooRealVar("mean2", "mean^2", 10, 0, 200)
51sigma = ROOT.RooRealVar("sigma", "sigma", 3, 0.1, 10)
52
53# Construct interpreted function mean = sqrt(mean^2)
54mean = ROOT.RooFormulaVar("mean", "mean", "sqrt(mean2)", [mean2])
55
56# Construct a gaussian g2(x,sqrt(mean2),sigma)
57g2 = ROOT.RooGaussian("g2", "h2", x, mean, sigma)
58
59# Generate toy data
60# ---------------------------------
61
62# Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian
63# dataset with mean 10 and width 3
64g1 = ROOT.RooGaussian("g1", "g1", x, ROOT.RooFit.RooConst(10), ROOT.RooFit.RooConst(3))
65data2 = g1.generate({x}, 1000)
66
67# Fit and plot tailored standard pdf
68# -------------------------------------------------------------------
69
70# Fit g2 to data from g1
71r = g2.fitTo(data2, Save=True) # ROOT.RooFitResult
72r.Print()
73
74# Plot data on frame and overlay projection of g2
75xframe2 = x.frame(Title="Tailored Gaussian pdf")
76data2.plotOn(xframe2)
77g2.plotOn(xframe2)
78
79# Draw all frames on a canvas
80c = ROOT.TCanvas("rf103_interprfuncs", "rf103_interprfuncs", 800, 400)
81c.Divide(2)
82c.cd(1)
83ROOT.gPad.SetLeftMargin(0.15)
84xframe.GetYaxis().SetTitleOffset(1.4)
85xframe.Draw()
86c.cd(2)
87ROOT.gPad.SetLeftMargin(0.15)
88xframe2.GetYaxis().SetTitleOffset(1.4)
89xframe2.Draw()
90
91c.SaveAs("rf103_interprfuncs.png")