Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf605_profilell.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605
6##
7## Working with the profile likelihood estimator
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 model and dataset
21# -----------------------------------------------
22
23# Observable
24x = ROOT.RooRealVar("x", "x", -20, 20)
25
26# Model (intentional strong correlations)
27mean = ROOT.RooRealVar("mean", "mean of g1 and g2", 0, -10, 10)
28sigma_g1 = ROOT.RooRealVar("sigma_g1", "width of g1", 3)
29g1 = ROOT.RooGaussian("g1", "g1", x, mean, sigma_g1)
30
31sigma_g2 = ROOT.RooRealVar("sigma_g2", "width of g2", 4, 3.0, 6.0)
32g2 = ROOT.RooGaussian("g2", "g2", x, mean, sigma_g2)
33
34frac = ROOT.RooRealVar("frac", "frac", 0.5, 0.0, 1.0)
35model = ROOT.RooAddPdf("model", "model", [g1, g2], [frac])
36
37# Generate 1000 events
38data = model.generate({x}, 1000)
39
40# Construct plain likelihood
41# ---------------------------------------------------
42
43# Construct unbinned likelihood
44nll = model.createNLL(data, NumCPU=2)
45
46# Minimize likelihood w.r.t all parameters before making plots
47ROOT.RooMinimizer(nll).migrad()
48
49# Plot likelihood scan frac
50frame1 = frac.frame(Bins=10, Range=(0.01, 0.95), Title="LL and profileLL in frac")
51nll.plotOn(frame1, ShiftToZero=True)
52
53# Plot likelihood scan in sigma_g2
54frame2 = sigma_g2.frame(Bins=10, Range=(3.3, 5.0), Title="LL and profileLL in sigma_g2")
55nll.plotOn(frame2, ShiftToZero=True)
56
57# Construct profile likelihood in frac
58# -----------------------------------------------------------------------
59
60# The profile likelihood estimator on nll for frac will minimize nll w.r.t
61# all floating parameters except frac for each evaluation
62
63pll_frac = nll.createProfile({frac})
64
65# Plot the profile likelihood in frac
66pll_frac.plotOn(frame1, LineColor="r")
67
68# Adjust frame maximum for visual clarity
69frame1.SetMinimum(0)
70frame1.SetMaximum(3)
71
72# Construct profile likelihood in sigma_g2
73# -------------------------------------------------------------------------------
74
75# The profile likelihood estimator on nll for sigma_g2 will minimize nll
76# w.r.t all floating parameters except sigma_g2 for each evaluation
77pll_sigmag2 = nll.createProfile({sigma_g2})
78
79# Plot the profile likelihood in sigma_g2
80pll_sigmag2.plotOn(frame2, LineColor="r")
81
82# Adjust frame maximum for visual clarity
83frame2.SetMinimum(0)
84frame2.SetMaximum(3)
85
86# Make canvas and draw ROOT.RooPlots
87c = ROOT.TCanvas("rf605_profilell", "rf605_profilell", 800, 400)
88c.Divide(2)
89c.cd(1)
90ROOT.gPad.SetLeftMargin(0.15)
91frame1.GetYaxis().SetTitleOffset(1.4)
92frame1.Draw()
93c.cd(2)
94ROOT.gPad.SetLeftMargin(0.15)
95frame2.GetYaxis().SetTitleOffset(1.4)
96frame2.Draw()
97
98c.SaveAs("rf605_profilell.png")