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