Logo ROOT   6.16/01
Reference Guide
rf607_fitresult.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
6##
7## Demonstration of options of the RooFitResult class
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15from __future__ import print_function
16import ROOT
17
18
19# Create pdf, data
20# --------------------------------
21
22# Declare observable x
23x = ROOT.RooRealVar("x", "x", 0, 10)
24
25# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
26# their parameters
27mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, -10, 10)
28sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5, 0.1, 10)
29sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1, 0.1, 10)
30
31sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
32sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
33
34# Build Chebychev polynomial p.d.f.
35a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
36a1 = ROOT.RooRealVar("a1", "a1", -0.2)
37bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
38
39# Sum the signal components into a composite signal p.d.f.
40sig1frac = ROOT.RooRealVar(
41 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
42sig = ROOT.RooAddPdf(
43 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
44
45# Sum the composite signal and background
46bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
47model = ROOT.RooAddPdf(
48 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
49
50# Generate 1000 events
51data = model.generate(ROOT.RooArgSet(x), 1000)
52
53# Fit pdf to data, save fit result
54# -------------------------------------------------------------
55
56# Perform fit and save result
57r = model.fitTo(data, ROOT.RooFit.Save())
58
59# Print fit results
60# ---------------------------------
61
62# Summary printing: Basic info plus final values of floating fit parameters
63r.Print()
64
65# Verbose printing: Basic info, of constant parameters, and
66# final values of floating parameters, correlations
67r.Print("v")
68
69# Visualize correlation matrix
70# -------------------------------------------------------
71
72# Construct 2D color plot of correlation matrix
73ROOT.gStyle.SetOptStat(0)
74ROOT.gStyle.SetPalette(1)
75hcorr = r.correlationHist()
76
77# Visualize ellipse corresponding to single correlation matrix element
78frame = ROOT.RooPlot(sigma1, sig1frac, 0.45, 0.60, 0.65, 0.90)
79frame.SetTitle("Covariance between sigma1 and sig1frac")
80r.plotOn(frame, sigma1, sig1frac, "ME12ABHV")
81
82# Access fit result information
83# ---------------------------------------------------------
84
85# Access basic information
86print("EDM = ", r.edm())
87print("-log(L) minimum = ", r.minNll())
88
89# Access list of final fit parameter values
90print("final value of floating parameters")
91r.floatParsFinal().Print("s")
92
93# Access correlation matrix elements
94print("correlation between sig1frac and a0 is ", r.correlation(
95 sig1frac, a0))
96print("correlation between bkgfrac and mean is ", r.correlation(
97 "bkgfrac", "mean"))
98
99# Extract covariance and correlation matrix as ROOT.TMatrixDSym
100cor = r.correlationMatrix()
101cov = r.covarianceMatrix()
102
103# Print correlation, matrix
104print("correlation matrix")
105cor.Print()
106print("covariance matrix")
107cov.Print()
108
109# Persist fit result in root file
110# -------------------------------------------------------------
111
112# Open ROOT file save save result
113f = ROOT.TFile("rf607_fitresult.root", "RECREATE")
114r.Write("rf607")
115f.Close()
116
117# In a clean ROOT session retrieve the persisted fit result as follows:
118# r = gDirectory.Get("rf607")
119
120c = ROOT.TCanvas("rf607_fitresult", "rf607_fitresult", 800, 400)
121c.Divide(2)
122c.cd(1)
123ROOT.gPad.SetLeftMargin(0.15)
124hcorr.GetYaxis().SetTitleOffset(1.4)
125hcorr.Draw("colz")
126c.cd(2)
127ROOT.gPad.SetLeftMargin(0.15)
128frame.GetYaxis().SetTitleOffset(1.6)
129frame.Draw()
130
131c.SaveAs("rf607_fitresult.png")
void Print(std::ostream &os, const OptionType &opt)