ROOT
Version v6.32
master
v6.34
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf601_intminuit.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
##
5
## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #601
6
##
7
## Interactive minimization with MINUIT
8
##
9
## \macro_image
10
## \macro_code
11
## \macro_output
12
##
13
## \date February 2018
14
## \authors Clemens Lange, Wouter Verkerke (C version)
15
16
17
import
ROOT
18
19
20
# Setup pdf and likelihood
21
# -----------------------------------------------
22
23
# Observable
24
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -20, 20)
25
26
# Model (intentional strong correlations)
27
mean =
ROOT.RooRealVar
(
"mean"
,
"mean of g1 and g2"
, 0)
28
sigma_g1 =
ROOT.RooRealVar
(
"sigma_g1"
,
"width of g1"
, 3)
29
g1 =
ROOT.RooGaussian
(
"g1"
,
"g1"
, x, mean, sigma_g1)
30
31
sigma_g2 =
ROOT.RooRealVar
(
"sigma_g2"
,
"width of g2"
, 4, 3.0, 6.0)
32
g2 =
ROOT.RooGaussian
(
"g2"
,
"g2"
, x, mean, sigma_g2)
33
34
frac =
ROOT.RooRealVar
(
"frac"
,
"frac"
, 0.5, 0.0, 1.0)
35
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [g1, g2], [frac])
36
37
# Generate 1000 events
38
data =
model.generate
({x}, 1000)
39
40
# Construct unbinned likelihood of model w.r.t. data
41
nll =
model.createNLL
(data)
42
43
# Interactive minimization, error analysis
44
# -------------------------------------------------------------------------------
45
46
# Create MINUIT interface object
47
m =
ROOT.RooMinimizer
(nll)
48
49
# Activate verbose logging of MINUIT parameter space stepping
50
m.setVerbose
(
True
)
51
52
# Call MIGRAD to minimize the likelihood
53
m.migrad
()
54
55
# Print values of all parameters, reflect values (and error estimates)
56
# that are back propagated from MINUIT
57
model.getParameters
({x}).
Print
(
"s"
)
58
59
# Disable verbose logging
60
m.setVerbose
(
False
)
61
62
# Run HESSE to calculate errors from d2L/dp2
63
m.hesse
()
64
65
# Print value (and error) of sigma_g2 parameter, reflects
66
# value and error back propagated from MINUIT
67
sigma_g2.Print
()
68
69
# Run MINOS on sigma_g2 parameter only
70
m.minos
({sigma_g2})
71
72
# Print value (and error) of sigma_g2 parameter, reflects
73
# value and error back propagated from MINUIT
74
sigma_g2.Print
()
75
76
# Saving results, contour plots
77
# ---------------------------------------------------------
78
79
# Save a snapshot of the fit result. ROOT.This object contains the initial
80
# fit parameters, final fit parameters, complete correlation
81
# matrix, EDM, minimized FCN , last MINUIT status code and
82
# the number of times the ROOT.RooFit function object has indicated evaluation
83
# problems (e.g. zero probabilities during likelihood evaluation)
84
r =
m.save
()
85
86
# Make contour plot of mx vs sx at 1,2, sigma
87
frame =
m.contour
(frac, sigma_g2, 1, 2, 3)
88
frame.SetTitle
(
"Contour plot"
)
89
90
# Print the fit result snapshot
91
r.Print
(
"v"
)
92
93
# Change parameter values, plotting
94
# -----------------------------------------------------------------
95
96
# At any moment you can manually change the value of a (constant)
97
# parameter
98
mean.setVal
(0.3)
99
100
# Rerun MIGRAD,HESSE
101
m.migrad
()
102
m.hesse
()
103
frac.Print
()
104
105
# Now fix sigma_g2
106
sigma_g2.setConstant
(
True
)
107
108
# Rerun MIGRAD,HESSE
109
m.migrad
()
110
m.hesse
()
111
frac.Print
()
112
113
c =
ROOT.TCanvas
(
"rf601_intminuit"
,
"rf601_intminuit"
, 600, 600)
114
ROOT.gPad.SetLeftMargin
(0.15)
115
frame.GetYaxis
().SetTitleOffset(1.4)
116
frame.Draw
()
117
118
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
rf601_intminuit.py
ROOT v6-32 - Reference Guide Generated on Wed Apr 2 2025 08:24:49 (GVA Time) using Doxygen 1.10.0