Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf705_linearmorph.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'SPECIAL PDFS' RooFit tutorial macro #705
6##
7## Linear interpolation between p.d.f shapes using the 'Alex Read' algorithm
8##
9## \macro_image
10## \macro_code
11## \macro_output
12##
13## \date February 2018
14## \authors Clemens Lange, Wouter Verkerke (C version)
15
16
17import ROOT
18
19
20# Create end point pdf shapes
21# ------------------------------------------------------
22
23# Observable
24x = ROOT.RooRealVar("x", "x", -20, 20)
25
26# Lower end point shape: a Gaussian
27g1mean = ROOT.RooRealVar("g1mean", "g1mean", -10)
28g1 = ROOT.RooGaussian("g1", "g1", x, g1mean, 2.0)
29
30# Upper end point shape: a Polynomial
31g2 = ROOT.RooPolynomial("g2", "g2", x, [-0.03, -0.001])
32
33# Create interpolating pdf
34# -----------------------------------------------
35
36# Create interpolation variable
37alpha = ROOT.RooRealVar("alpha", "alpha", 0, 1.0)
38
39# Specify sampling density on observable and interpolation variable
40x.setBins(1000, "cache")
41alpha.setBins(50, "cache")
42
43# Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
44# and g2(x) at a=a_max
45lmorph = ROOT.RooIntegralMorph("lmorph", "lmorph", g1, g2, x, alpha)
46
47# Plot interpolating pdf aat various alphas a l p h a
48# -----------------------------------------------------------------------------
49
50# Show end points as blue curves
51frame1 = x.frame()
52g1.plotOn(frame1)
53g2.plotOn(frame1)
54
55# Show interpolated shapes in red
56alpha.setVal(0.125)
57lmorph.plotOn(frame1, LineColor="r")
58alpha.setVal(0.25)
59lmorph.plotOn(frame1, LineColor="r")
60alpha.setVal(0.375)
61lmorph.plotOn(frame1, LineColor="r")
62alpha.setVal(0.50)
63lmorph.plotOn(frame1, LineColor="r")
64alpha.setVal(0.625)
65lmorph.plotOn(frame1, LineColor="r")
66alpha.setVal(0.75)
67lmorph.plotOn(frame1, LineColor="r")
68alpha.setVal(0.875)
69lmorph.plotOn(frame1, LineColor="r")
70alpha.setVal(0.95)
71lmorph.plotOn(frame1, LineColor="r")
72
73# Show 2D distribution of pdf(x,alpha)
74# -----------------------------------------------------------------------
75
76# Create 2D histogram
77hh = lmorph.createHistogram("hh", x, Binning=40, YVar=dict(var=alpha, Binning=40))
78hh.SetLineColor(ROOT.kBlue)
79
80# Fit pdf to dataset with alpha=0.8
81# -----------------------------------------------------------------
82
83# Generate a toy dataset alpha = 0.8
84alpha.setVal(0.8)
85data = lmorph.generate({x}, 1000)
86
87# Fit pdf to toy data
88lmorph.setCacheAlpha(True)
89lmorph.fitTo(data, Verbose=True, PrintLevel=-1)
90
91# Plot fitted pdf and data overlaid
92frame2 = x.frame(Bins=100)
93data.plotOn(frame2)
94lmorph.plotOn(frame2)
95
96# Scan -log(L) vs alpha
97# -----------------------------------------
98
99# Show scan -log(L) of dataset w.r.t alpha
100frame3 = alpha.frame(Bins=100, Range=(0.1, 0.9))
101
102# Make 2D pdf of histogram
103nll = lmorph.createNLL(data)
104nll.plotOn(frame3, ShiftToZero=True)
105
106lmorph.setCacheAlpha(False)
107
108c = ROOT.TCanvas("rf705_linearmorph", "rf705_linearmorph", 800, 800)
109c.Divide(2, 2)
110c.cd(1)
111ROOT.gPad.SetLeftMargin(0.15)
112frame1.GetYaxis().SetTitleOffset(1.6)
113frame1.Draw()
114c.cd(2)
115ROOT.gPad.SetLeftMargin(0.20)
116hh.GetZaxis().SetTitleOffset(2.5)
117hh.Draw("surf")
118c.cd(3)
119ROOT.gPad.SetLeftMargin(0.15)
120frame3.GetYaxis().SetTitleOffset(1.4)
121frame3.Draw()
122c.cd(4)
123ROOT.gPad.SetLeftMargin(0.15)
124frame2.GetYaxis().SetTitleOffset(1.4)
125frame2.Draw()
126
127c.SaveAs("rf705_linearmorph.png")