Rf 6 0 5_Profilell

'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605

Working with the profile likelihood estimator

Author: Wouter Verkerke (C version)
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Monday, January 17, 2022 at 10:08 AM.

In [ ]:
import ROOT

Create model and dataset

Observable

In [ ]:
x = ROOT.RooRealVar("x", "x", -20, 20)

Model (intentional strong correlations)

In [ ]:
mean = ROOT.RooRealVar("mean", "mean of g1 and g2", 0, -10, 10)
sigma_g1 = ROOT.RooRealVar("sigma_g1", "width of g1", 3)
g1 = ROOT.RooGaussian("g1", "g1", x, mean, sigma_g1)

sigma_g2 = ROOT.RooRealVar("sigma_g2", "width of g2", 4, 3.0, 6.0)
g2 = ROOT.RooGaussian("g2", "g2", x, mean, sigma_g2)

frac = ROOT.RooRealVar("frac", "frac", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "model", [g1, g2], [frac])

Generate 1000 events

In [ ]:
data = model.generate({x}, 1000)

Construct plain likelihood

Construct unbinned likelihood

In [ ]:
nll = model.createNLL(data, NumCPU=2)

Minimize likelihood w.r.t all parameters before making plots

In [ ]:
ROOT.RooMinimizer(nll).migrad()

Plot likelihood scan frac

In [ ]:
frame1 = frac.frame(Bins=10, Range=(0.01, 0.95), Title="LL and profileLL in frac")
nll.plotOn(frame1, ShiftToZero=True)

Plot likelihood scan in sigma_g2

In [ ]:
frame2 = sigma_g2.frame(Bins=10, Range=(3.3, 5.0), Title="LL and profileLL in sigma_g2")
nll.plotOn(frame2, ShiftToZero=True)

Construct profile likelihood in frac

The profile likelihood estimator on nll for frac will minimize nll w.r.t all floating parameters except frac for each evaluation

In [ ]:
pll_frac = nll.createProfile({frac})

Plot the profile likelihood in frac

In [ ]:
pll_frac.plotOn(frame1, LineColor="r")

Adjust frame maximum for visual clarity

In [ ]:
frame1.SetMinimum(0)
frame1.SetMaximum(3)

Construct profile likelihood in sigma_g2

The profile likelihood estimator on nll for sigma_g2 will minimize nll w.r.t all floating parameters except sigma_g2 for each evaluation

In [ ]:
pll_sigmag2 = nll.createProfile({sigma_g2})

Plot the profile likelihood in sigma_g2

In [ ]:
pll_sigmag2.plotOn(frame2, LineColor="r")

Adjust frame maximum for visual clarity

In [ ]:
frame2.SetMinimum(0)
frame2.SetMaximum(3)

Make canvas and draw ROOT.RooPlots

In [ ]:
c = ROOT.TCanvas("rf605_profilell", "rf605_profilell", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.4)
frame1.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.4)
frame2.Draw()

c.SaveAs("rf605_profilell.png")

Draw all canvases

In [ ]:
from ROOT import gROOT 
gROOT.GetListOfCanvases().Draw()