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