Logo ROOT   6.16/01
Reference Guide
rf706_histpdf.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'SPECIAL PDFS' RooFit tutorial macro #706
6##
7## Histogram based p.d.f.s and functions
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13
14
15import ROOT
16
17
18# Create pdf for sampling
19# ---------------------------------------------
20
21x = ROOT.RooRealVar("x", "x", 0, 20)
22p = ROOT.RooPolynomial("p", "p", x, ROOT.RooArgList(ROOT.RooFit.RooConst(
23 0.01), ROOT.RooFit.RooConst(-0.01), ROOT.RooFit.RooConst(0.0004)))
24
25# Create low stats histogram
26# ---------------------------------------------------
27
28# Sample 500 events from p
29x.setBins(20)
30data1 = p.generate(ROOT.RooArgSet(x), 500)
31
32# Create a binned dataset with 20 bins and 500 events
33hist1 = data1.binnedClone()
34
35# Represent data in dh as pdf in x
36histpdf1 = ROOT.RooHistPdf("histpdf1", "histpdf1", ROOT.RooArgSet(x), hist1, 0)
37
38# Plot unbinned data and histogram pdf overlaid
39frame1 = x.frame(ROOT.RooFit.Title(
40 "Low statistics histogram pdf"), ROOT.RooFit.Bins(100))
41data1.plotOn(frame1)
42histpdf1.plotOn(frame1)
43
44# Create high stats histogram
45# -----------------------------------------------------
46
47# Sample 100000 events from p
48x.setBins(10)
49data2 = p.generate(ROOT.RooArgSet(x), 100000)
50
51# Create a binned dataset with 10 bins and 100K events
52hist2 = data2.binnedClone()
53
54# Represent data in dh as pdf in x, 2nd order interpolation
55histpdf2 = ROOT.RooHistPdf("histpdf2", "histpdf2", ROOT.RooArgSet(x), hist2, 2)
56
57# Plot unbinned data and histogram pdf overlaid
58frame2 = x.frame(ROOT.RooFit.Title(
59 "High stats histogram pdf with interpolation"), ROOT.RooFit.Bins(100))
60data2.plotOn(frame2)
61histpdf2.plotOn(frame2)
62
63c = ROOT.TCanvas("rf706_histpdf", "rf706_histpdf", 800, 400)
64c.Divide(2)
65c.cd(1)
66ROOT.gPad.SetLeftMargin(0.15)
67frame1.GetYaxis().SetTitleOffset(1.4)
68frame1.Draw()
69c.cd(2)
70ROOT.gPad.SetLeftMargin(0.15)
71frame2.GetYaxis().SetTitleOffset(1.8)
72frame2.Draw()
73
74c.SaveAs("rf706_histpdf.png")