Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf601_intminuit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #601
5##
6## Interactive minimization with MINUIT
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# Setup pdf and likelihood
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)
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 unbinned likelihood of model w.r.t. data
40nll = model.createNLL(data)
41
42# Interactive minimization, error analysis
43# -------------------------------------------------------------------------------
44
45# Create MINUIT interface object
46m = ROOT.RooMinimizer(nll)
47
48# Activate verbose logging of MINUIT parameter space stepping
49m.setVerbose(True)
50
51# Call MIGRAD to minimize the likelihood
52m.migrad()
53
54# Print values of all parameters, reflect values (and error estimates)
55# that are back propagated from MINUIT
56model.getParameters({x}).Print("s")
57
58# Disable verbose logging
59m.setVerbose(False)
60
61# Run HESSE to calculate errors from d2L/dp2
62m.hesse()
63
64# Print value (and error) of sigma_g2 parameter, reflects
65# value and error back propagated from MINUIT
66sigma_g2.Print()
67
68# Run MINOS on sigma_g2 parameter only
69m.minos({sigma_g2})
70
71# Print value (and error) of sigma_g2 parameter, reflects
72# value and error back propagated from MINUIT
73sigma_g2.Print()
74
75# Saving results, contour plots
76# ---------------------------------------------------------
77
78# Save a snapshot of the fit result. ROOT.This object contains the initial
79# fit parameters, final fit parameters, complete correlation
80# matrix, EDM, minimized FCN , last MINUIT status code and
81# the number of times the ROOT.RooFit function object has indicated evaluation
82# problems (e.g. zero probabilities during likelihood evaluation)
83r = m.save()
84
85# Make contour plot of mx vs sx at 1,2, sigma
86frame = m.contour(frac, sigma_g2, 1, 2, 3)
87frame.SetTitle("Contour plot")
88
89# Print the fit result snapshot
90r.Print("v")
91
92# Change parameter values, plotting
93# -----------------------------------------------------------------
94
95# At any moment you can manually change the value of a (constant)
96# parameter
97mean.setVal(0.3)
98
99# Rerun MIGRAD,HESSE
100m.migrad()
101m.hesse()
102frac.Print()
103
104# Now fix sigma_g2
105sigma_g2.setConstant(True)
106
107# Rerun MIGRAD,HESSE
108m.migrad()
109m.hesse()
110frac.Print()
111
112c = ROOT.TCanvas("rf601_intminuit", "rf601_intminuit", 600, 600)
113ROOT.gPad.SetLeftMargin(0.15)
114frame.GetYaxis().SetTitleOffset(1.4)
115frame.Draw()
116
117c.SaveAs("rf601_intminuit.png")
void Print(GNN_Data &d, std::string txt="")