Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf705_linearmorph.py File Reference

Namespaces

namespace  rf705_linearmorph
 

Detailed Description

View in nbviewer Open in SWAN

'SPECIAL PDFS' RooFit tutorial macro #705

Linear interpolation between p.d.f shapes using the 'Alex Read' algorithm

import ROOT
# Create end point pdf shapes
# ------------------------------------------------------
# Observable
x = ROOT.RooRealVar("x", "x", -20, 20)
# Lower end point shape: a Gaussian
g1mean = ROOT.RooRealVar("g1mean", "g1mean", -10)
g1 = ROOT.RooGaussian("g1", "g1", x, g1mean, ROOT.RooFit.RooConst(2))
# Upper end point shape: a Polynomial
g2 = ROOT.RooPolynomial("g2", "g2", x, [-0.03, -0.001])
# Create interpolating pdf
# -----------------------------------------------
# Create interpolation variable
alpha = ROOT.RooRealVar("alpha", "alpha", 0, 1.0)
# Specify sampling density on observable and interpolation variable
x.setBins(1000, "cache")
alpha.setBins(50, "cache")
# Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
# and g2(x) at a=a_max
lmorph = ROOT.RooIntegralMorph("lmorph", "lmorph", g1, g2, x, alpha)
# Plot interpolating pdf aat various alphas a l p h a
# -----------------------------------------------------------------------------
# Show end points as blue curves
frame1 = x.frame()
g1.plotOn(frame1)
g2.plotOn(frame1)
# Show interpolated shapes in red
alpha.setVal(0.125)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.25)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.375)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.50)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.625)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.75)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.875)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.95)
lmorph.plotOn(frame1, LineColor="r")
# Show 2D distribution of pdf(x,alpha)
# -----------------------------------------------------------------------
# Create 2D histogram
hh = lmorph.createHistogram("hh", x, Binning=40, YVar=dict(var=alpha, Binning=40))
hh.SetLineColor(ROOT.kBlue)
# Fit pdf to dataset with alpha=0.8
# -----------------------------------------------------------------
# Generate a toy dataset alpha = 0.8
alpha.setVal(0.8)
data = lmorph.generate({x}, 1000)
# Fit pdf to toy data
lmorph.setCacheAlpha(True)
lmorph.fitTo(data, Verbose=True, PrintLevel=-1)
# Plot fitted pdf and data overlaid
frame2 = x.frame(Bins=100)
data.plotOn(frame2)
lmorph.plotOn(frame2)
# Scan -log(L) vs alpha
# -----------------------------------------
# Show scan -log(L) of dataset w.r.t alpha
frame3 = alpha.frame(Bins=100, Range=(0.1, 0.9))
# Make 2D pdf of histogram
nll = lmorph.createNLL(data)
nll.plotOn(frame3, ShiftToZero=True)
lmorph.setCacheAlpha(False)
c = ROOT.TCanvas("rf705_linearmorph", "rf705_linearmorph", 800, 800)
c.Divide(2, 2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.6)
frame1.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.20)
hh.GetZaxis().SetTitleOffset(2.5)
hh.Draw("surf")
c.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
frame3.GetYaxis().SetTitleOffset(1.4)
frame3.Draw()
c.cd(4)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.4)
frame2.Draw()
c.SaveAs("rf705_linearmorph.png")
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x67fd470 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x6840e20 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x_alpha for nset (x,alpha) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x64d5910 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x68ccc50 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0 from preexisting content.
[#0] PROGRESS:Eval -- RooIntegralMorph::fillCacheObject(lmorph) filling multi-dimensional cache..................................................
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x27485b0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (g1,g2)
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for alpha: using 0.1
prevFCN = 2639.749221 alpha=0.8026,
prevFCN = 2639.409749 alpha=0.7974,
prevFCN = 2640.17258 alpha=0.8003,
prevFCN = 2639.703144 alpha=0.7997,
prevFCN = 2639.796557 alpha=0.8116,
prevFCN = 2638.919064 alpha=0.8119,
prevFCN = 2638.930985 alpha=0.8113,
prevFCN = 2638.908768 alpha=0.8094,
prevFCN = 2638.901287 alpha=0.8102,
prevFCN = 2638.885143 alpha=0.8103,
prevFCN = 2638.886672 alpha=0.8101,
prevFCN = 2638.883909 alpha=0.8095,
prevFCN = 2638.898442 alpha=0.8098,
prevFCN = 2638.889204 alpha=0.8099,
prevFCN = 2638.884962 alpha=0.8104,
prevFCN = 2638.887231 alpha=0.8098,
prevFCN = 2638.888592 alpha=0.8102,
prevFCN = 2638.885076 alpha=0.81,
prevFCN = 2638.883889 alpha=0.8098,
prevFCN = 2638.888745 alpha=0.81,
prevFCN = 2638.883589 alpha=0.81,
prevFCN = 2638.883345 alpha=0.8095,
prevFCN = 2638.898453 alpha=0.8098,
prevFCN = 2638.89011 alpha=0.8099,
prevFCN = 2638.886243 alpha=0.81,
prevFCN = 2638.884385 alpha=0.81,
prevFCN = 2638.883475 alpha=0.8101,
prevFCN = 2638.884457 alpha=0.8099,
prevFCN = 2638.885715 alpha=0.8101,
prevFCN = 2638.883988 alpha=0.81,
prevFCN = 2638.884446 alpha=0.81,
prevFCN = 2638.883552 alpha=0.81,
prevFCN = 2638.883372 alpha=0.81,
prevFCN = 2638.883349 alpha=0.81,
prevFCN = 2638.883345 alpha=0.81,
prevFCN = 2638.883345 alpha=0.8101,
prevFCN = 2638.883988 alpha=0.81,
prevFCN = 2638.884446 alpha=0.81,
prevFCN = 2638.883467 alpha=0.81,
prevFCN = 2638.883226 alpha=0.81,
prevFCN = 2638.883369 alpha=0.81,
prevFCN = 2638.883321 alpha=0.81,
prevFCN = 2638.883345 alpha=0.81,
prevFCN = 2638.883369 alpha=0.81,
prevFCN = 2638.883321 alpha=0.8103,
prevFCN = 2638.886457 alpha=0.8097,
prevFCN = 2638.89054 alpha=0.8131,
prevFCN = 2638.992315 alpha=0.8069,
prevFCN = 2639.023618 alpha=0.8171,
prevFCN = 2639.367946 alpha=0.8029,
prevFCN = 2639.374671 alpha=0.8102,
prevFCN = 2638.885206 alpha=0.8098,
prevFCN = 2638.887615 alpha=0.8101,
prevFCN = 2638.883674 alpha=0.81,
prevFCN = 2638.883555 alpha=0.81, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x2743110 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x6fae720 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x5ea8680 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x5ea8680 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C version)

Definition in file rf705_linearmorph.py.