Likelihood and minimization: demonstration of options of the RooFitResult class 
 
  
from __future__ import print_function
import ROOT
 
 
 
x = ROOT.RooRealVar("x", "x", 0, 10)
 
mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, -10, 10)
sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5, 0.1, 10)
sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1, 0.1, 10)
 
sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
 
a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
a1 = ROOT.RooRealVar("a1", "a1", -0.2)
bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
 
sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
 
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
 
data = model.generate({x}, 1000)
 
 
r = model.fitTo(data, Save=True, PrintLevel=-1)
 
 
r.Print()
 
r.Print("v")
 
 
ROOT.gStyle.SetOptStat(0)
ROOT.gStyle.SetPalette(1)
hcorr = r.correlationHist()
 
frame = ROOT.RooPlot(sigma1, sig1frac, 0.45, 0.60, 0.65, 0.90)
frame.SetTitle("Covariance between sigma1 and sig1frac")
r.plotOn(frame, sigma1, sig1frac, "ME12ABHV")
 
 
print("EDM = ", r.edm())
print("-log(L) minimum = ", r.minNll())
 
print("final value of floating parameters")
r.floatParsFinal().Print("s")
 
print("correlation between sig1frac and a0 is  ", r.correlation(sig1frac, a0))
print("correlation between bkgfrac and mean is ", r.correlation("bkgfrac", "mean"))
 
cor = r.correlationMatrix()
cov = r.covarianceMatrix()
 
print("correlation matrix")
cor.Print()
print("covariance matrix")
cov.Print()
 
 
f = ROOT.TFile("rf607_fitresult.root", "RECREATE")
r.Write("rf607")
f.Close()
 
 
c = ROOT.TCanvas("rf607_fitresult", "rf607_fitresult", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
hcorr.GetYaxis().SetTitleOffset(1.4)
hcorr.Draw("colz")
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.6)
frame.Draw()
 
c.SaveAs("rf607_fitresult.png")
  [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --  The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
 
  RooFitResult: minimized FCN value: 1885.34, estimated distance to minimum: 0.000205499
                covariance matrix quality: Full, accurate covariance matrix
                Status : MINIMIZE=0 HESSE=0 
 
    Floating Parameter    FinalValue +/-  Error   
  --------------------  --------------------------
                    a0    7.2825e-01 +/-  1.11e-01
               bkgfrac    4.3439e-01 +/-  8.36e-02
                  mean    5.0346e+00 +/-  3.36e-02
              sig1frac    7.7835e-01 +/-  9.70e-02
                sigma1    5.2340e-01 +/-  4.51e-02
                sigma2    1.7767e+00 +/-  1.16e+00
 
 
  RooFitResult: minimized FCN value: 1885.34, estimated distance to minimum: 0.000205499
                covariance matrix quality: Full, accurate covariance matrix
                Status : MINIMIZE=0 HESSE=0 
 
    Constant Parameter    Value     
  --------------------  ------------
                    a1   -2.0000e-01
 
    Floating Parameter  InitialValue    FinalValue +/-  Error     GblCorr.
  --------------------  ------------  --------------------------  --------
                    a0    5.0000e-01    7.2825e-01 +/-  1.11e-01  <none>
               bkgfrac    5.0000e-01    4.3439e-01 +/-  8.36e-02  <none>
                  mean    5.0000e+00    5.0346e+00 +/-  3.36e-02  <none>
              sig1frac    8.0000e-01    7.7835e-01 +/-  9.70e-02  <none>
                sigma1    5.0000e-01    5.2340e-01 +/-  4.51e-02  <none>
                sigma2    1.0000e+00    1.7767e+00 +/-  1.16e+00  <none>
 
  1) RooRealVar::       a0 = 0.728245 +/- 0.111109
  2) RooRealVar::  bkgfrac = 0.434386 +/- 0.0836079
  3) RooRealVar::     mean = 5.03463 +/- 0.0336219
  4) RooRealVar:: sig1frac = 0.778347 +/- 0.0969912
  5) RooRealVar::   sigma1 = 0.523396 +/- 0.0451307
  6) RooRealVar::   sigma2 = 1.77668 +/- 1.15533
 
6x6 matrix is as follows
 
     |      0    |      1    |      2    |      3    |      4    |
----------------------------------------------------------------------
   0 |          1     -0.7952    -0.02552     -0.3779      0.4111 
   1 |    -0.7952           1    -0.05102      0.6023     -0.3876 
   2 |   -0.02552    -0.05102           1     -0.0873    -0.04206 
   3 |    -0.3779      0.6023     -0.0873           1      0.2966 
   4 |     0.4111     -0.3876    -0.04206      0.2966           1 
   5 |     0.8272     -0.8708     0.01245     -0.2609      0.5799 
 
 
     |      5    |
----------------------------------------------------------------------
   0 |     0.8272 
   1 |    -0.8708 
   2 |    0.01245 
   3 |    -0.2609 
   4 |     0.5799 
   5 |          1 
 
 
6x6 matrix is as follows
 
     |      0    |      1    |      2    |      3    |      4    |
----------------------------------------------------------------------
   0 |    0.01261   -0.007502  -9.635e-05   -0.004154    0.002084 
   1 |  -0.007502    0.007058  -0.0001441    0.004954    -0.00147 
   2 | -9.635e-05  -0.0001441     0.00113  -0.0002873  -6.383e-05 
   3 |  -0.004154    0.004954  -0.0002873    0.009583     0.00131 
   4 |   0.002084    -0.00147  -6.383e-05     0.00131    0.002037 
   5 |     0.1091    -0.08595   0.0004916       -0.03     0.03075 
 
 
     |      5    |
----------------------------------------------------------------------
   0 |     0.1091 
   1 |   -0.08595 
   2 |  0.0004916 
   3 |      -0.03 
   4 |    0.03075 
   5 |       1.38 
 
EDM =  0.00020549915143065197
-log(L) minimum =  1885.343815305069
final value of floating parameters
correlation between sig1frac and a0 is   -0.3778511951671703
correlation between bkgfrac and mean is  -0.05102316430532121
correlation matrix
covariance matrix
- Date
 - February 2018 
 
- Authors
 - Clemens Lange, Wouter Verkerke (C++ version) 
 
Definition in file rf607_fitresult.py.