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