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., 1.)
32a1 = ROOT.RooRealVar("a1", "a1", -0.2)
33bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
34
35# Sum the signal components into a composite signal pdf
36sig1frac = ROOT.RooRealVar(
37 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
38sig = ROOT.RooAddPdf(
39 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
40
41# Sum the composite signal and background
42bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
43model = ROOT.RooAddPdf(
44 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
45
46# Generate 1000 events
47data = model.generate(ROOT.RooArgSet(x), 1000)
48
49# Fit pdf to data, save fit result
50# -------------------------------------------------------------
51
52# Perform fit and save result
53r = model.fitTo(data, ROOT.RooFit.Save())
54
55# Print fit results
56# ---------------------------------
57
58# Summary printing: Basic info plus final values of floating fit parameters
59r.Print()
60
61# Verbose printing: Basic info, of constant parameters, and
62# final values of floating parameters, correlations
63r.Print("v")
64
65# Visualize correlation matrix
66# -------------------------------------------------------
67
68# Construct 2D color plot of correlation matrix
69ROOT.gStyle.SetOptStat(0)
70ROOT.gStyle.SetPalette(1)
71hcorr = r.correlationHist()
72
73# Visualize ellipse corresponding to single correlation matrix element
74frame = ROOT.RooPlot(sigma1, sig1frac, 0.45, 0.60, 0.65, 0.90)
75frame.SetTitle("Covariance between sigma1 and sig1frac")
76r.plotOn(frame, sigma1, sig1frac, "ME12ABHV")
77
78# Access fit result information
79# ---------------------------------------------------------
80
81# Access basic information
82print("EDM = ", r.edm())
83print("-log(L) minimum = ", r.minNll())
84
85# Access list of final fit parameter values
86print("final value of floating parameters")
87r.floatParsFinal().Print("s")
88
89# Access correlation matrix elements
90print("correlation between sig1frac and a0 is ", r.correlation(
91 sig1frac, a0))
92print("correlation between bkgfrac and mean is ", r.correlation(
93 "bkgfrac", "mean"))
94
95# Extract covariance and correlation matrix as ROOT.TMatrixDSym
96cor = r.correlationMatrix()
97cov = r.covarianceMatrix()
98
99# Print correlation, matrix
100print("correlation matrix")
101cor.Print()
102print("covariance matrix")
103cov.Print()
104
105# Persist fit result in root file
106# -------------------------------------------------------------
107
108# Open ROOT file save save result
109f = ROOT.TFile("rf607_fitresult.root", "RECREATE")
110r.Write("rf607")
111f.Close()
112
113# In a clean ROOT session retrieve the persisted fit result as follows:
114# r = gDirectory.Get("rf607")
115
116c = ROOT.TCanvas("rf607_fitresult", "rf607_fitresult", 800, 400)
117c.Divide(2)
118c.cd(1)
119ROOT.gPad.SetLeftMargin(0.15)
120hcorr.GetYaxis().SetTitleOffset(1.4)
121hcorr.Draw("colz")
122c.cd(2)
123ROOT.gPad.SetLeftMargin(0.15)
124frame.GetYaxis().SetTitleOffset(1.6)
125frame.Draw()
126
127c.SaveAs("rf607_fitresult.png")