ROOT
master
Reference Guide
Loading...
Searching...
No Matches
rf601_intminuit.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
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
16
import
ROOT
17
18
19
# Setup pdf and likelihood
20
# -----------------------------------------------
21
22
# Observable
23
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -20, 20)
24
25
# Model (intentional strong correlations)
26
mean =
ROOT.RooRealVar
(
"mean"
,
"mean of g1 and g2"
, 0)
27
sigma_g1 =
ROOT.RooRealVar
(
"sigma_g1"
,
"width of g1"
, 3)
28
g1 =
ROOT.RooGaussian
(
"g1"
,
"g1"
, x, mean, sigma_g1)
29
30
sigma_g2 =
ROOT.RooRealVar
(
"sigma_g2"
,
"width of g2"
, 4, 3.0, 6.0)
31
g2 =
ROOT.RooGaussian
(
"g2"
,
"g2"
, x, mean, sigma_g2)
32
33
frac =
ROOT.RooRealVar
(
"frac"
,
"frac"
, 0.5, 0.0, 1.0)
34
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [g1, g2], [frac])
35
36
# Generate 1000 events
37
data =
model.generate
({x}, 1000)
38
39
# Construct unbinned likelihood of model w.r.t. data
40
nll =
model.createNLL
(data)
41
42
# Interactive minimization, error analysis
43
# -------------------------------------------------------------------------------
44
45
# Create MINUIT interface object
46
m =
ROOT.RooMinimizer
(nll)
47
48
# Activate verbose logging of MINUIT parameter space stepping
49
m.setVerbose
(
True
)
50
51
# Call MIGRAD to minimize the likelihood
52
m.migrad
()
53
54
# Print values of all parameters, reflect values (and error estimates)
55
# that are back propagated from MINUIT
56
model.getParameters
({x}).
Print
(
"s"
)
57
58
# Disable verbose logging
59
m.setVerbose
(
False
)
60
61
# Run HESSE to calculate errors from d2L/dp2
62
m.hesse
()
63
64
# Print value (and error) of sigma_g2 parameter, reflects
65
# value and error back propagated from MINUIT
66
sigma_g2.Print
()
67
68
# Run MINOS on sigma_g2 parameter only
69
m.minos
({sigma_g2})
70
71
# Print value (and error) of sigma_g2 parameter, reflects
72
# value and error back propagated from MINUIT
73
sigma_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)
83
r =
m.save
()
84
85
# Make contour plot of mx vs sx at 1,2, sigma
86
frame =
m.contour
(frac, sigma_g2, 1, 2, 3)
87
frame.SetTitle
(
"Contour plot"
)
88
89
# Print the fit result snapshot
90
r.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
97
mean.setVal
(0.3)
98
99
# Rerun MIGRAD,HESSE
100
m.migrad
()
101
m.hesse
()
102
frac.Print
()
103
104
# Now fix sigma_g2
105
sigma_g2.setConstant
(
True
)
106
107
# Rerun MIGRAD,HESSE
108
m.migrad
()
109
m.hesse
()
110
frac.Print
()
111
112
c =
ROOT.TCanvas
(
"rf601_intminuit"
,
"rf601_intminuit"
, 600, 600)
113
ROOT.gPad.SetLeftMargin
(0.15)
114
frame.GetYaxis
().SetTitleOffset(1.4)
115
frame.Draw
()
116
117
c.SaveAs
(
"rf601_intminuit.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
Print
void Print(GNN_Data &d, std::string txt="")
Definition
TMVA_SOFIE_GNN_Application.C:59
tutorials
roofit
roofit
rf601_intminuit.py
ROOT master - Reference Guide Generated on Thu Mar 27 2025 04:32:53 (GVA Time) using Doxygen 1.10.0