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