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