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