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