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

Namespaces

namespace  rf605_profilell
 

Detailed Description

View in nbviewer Open in SWAN

'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605

Working with the profile likelihood estimator

import ROOT
# Create model and dataset
# -----------------------------------------------
# Observable
x = ROOT.RooRealVar("x", "x", -20, 20)
# Model (intentional strong correlations)
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
data = model.generate({x}, 1000)
# Construct plain likelihood
# ---------------------------------------------------
# Construct unbinned likelihood
nll = model.createNLL(data, NumCPU=2)
# Minimize likelihood w.r.t all parameters before making plots
ROOT.RooMinimizer(nll).migrad()
# Plot likelihood scan frac
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
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
pll_frac = nll.createProfile({frac})
# Plot the profile likelihood in frac
pll_frac.plotOn(frame1, LineColor="r")
# Adjust frame maximum for visual clarity
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
pll_sigmag2 = nll.createProfile({sigma_g2})
# Plot the profile likelihood in sigma_g2
pll_sigmag2.plotOn(frame2, LineColor="r")
# Adjust frame maximum for visual clarity
frame2.SetMinimum(0)
frame2.SetMaximum(3)
# Make canvas and draw ROOT.RooPlots
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")
Int_t migrad()
[#0] WARNING:InputArguments -- The parameter 'sigma_g1' with range [-1e+30, 1e+30] of the RooGaussian 'g1' exceeds the safe range of (0, inf). Advise to limit its range.
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 frac 5.00000e-01 1.00000e-01 0.00000e+00 1.00000e+00
2 mean 0.00000e+00 2.00000e+00 -1.00000e+01 1.00000e+01
3 sigma_g2 4.00000e+00 3.00000e-01 3.00000e+00 6.00000e+00
**********
** 3 **SET ERR 0.5
**********
**********
** 4 **SET PRINT 1
**********
**********
** 5 **SET STR 1
**********
NOW USING STRATEGY 1: TRY TO BALANCE SPEED AGAINST RELIABILITY
**********
** 6 **MIGRAD 1500 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 2 remote server process.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=2660.22 FROM MIGRAD STATUS=INITIATE 10 CALLS 11 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 frac 5.00000e-01 1.00000e-01 2.01358e-01 -5.61980e+00
2 mean 0.00000e+00 2.00000e+00 2.01358e-01 -7.16779e+00
3 sigma_g2 4.00000e+00 3.00000e-01 2.14402e-01 7.28535e+00
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=2659.74 FROM MIGRAD STATUS=CONVERGED 67 CALLS 68 TOTAL
EDM=5.19798e-06 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 frac 6.23972e-01 1.64510e-01 5.33134e-03 6.83204e-03
2 mean 4.57491e-03 1.09369e-01 3.87767e-04 -1.84350e-01
3 sigma_g2 4.11576e+00 4.07375e-01 4.33560e-03 -6.97269e-03
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 3 ERR DEF=0.5
2.817e-02 -1.610e-03 6.258e-02
-1.610e-03 1.196e-02 -4.302e-03
6.258e-02 -4.302e-03 1.705e-01
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3
1 0.90293 1.000 -0.088 0.903
2 0.09533 -0.088 1.000 -0.095
3 0.90308 0.903 -0.095 1.000
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 2 remote server process.
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 2 remote server process.
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[frac]) Creating instance of MINUIT
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[frac]) determining minimum likelihood for current configurations w.r.t all observable
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 2 remote server process.
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[frac]) minimum found at (frac=0.623915)
..................................................................................
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[sigma_g2]) Creating instance of MINUIT
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[sigma_g2]) determining minimum likelihood for current configurations w.r.t all observable
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 2 remote server process.
[#1] INFO:Minimization -- RooProfileLL::evaluate(nll_model_modelData_Profile[sigma_g2]) minimum found at (sigma_g2=4.11588)
....................................................................................
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C version)

Definition in file rf605_profilell.py.