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