ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
// 
// Demonstration of options of the RooFitResult class
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooAddPdf.h"
#include "RooChebychev.h"
#include "RooFitResult.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "TFile.h"
#include "TStyle.h"
#include "TH2.h"
#include "TMatrixDSym.h"

using namespace RooFit ;


void rf607_fitresult()
{
  // C r e a t e   p d f ,   d a t a
  // --------------------------------

  // Declare observable x
  RooRealVar x("x","x",0,10) ;

  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their paramaters
  RooRealVar mean("mean","mean of gaussians",5,-10,10) ;
  RooRealVar sigma1("sigma1","width of gaussians",0.5,0.1,10) ;
  RooRealVar sigma2("sigma2","width of gaussians",1,0.1,10) ;

  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;  
  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;  
  
  // Build Chebychev polynomial p.d.f.  
  RooRealVar a0("a0","a0",0.5,0.,1.) ;
  RooRealVar a1("a1","a1",-0.2) ;
  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;

  // Sum the signal components into a composite signal p.d.f.
  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;

  // Sum the composite signal and background 
  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
  RooAddPdf  model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;

  // Generate 1000 events
  RooDataSet* data = model.generate(x,1000) ;



  // F i t   p d f   t o   d a t a ,   s a v e   f i t r e s u l t 
  // -------------------------------------------------------------

  // Perform fit and save result
  RooFitResult* r = model.fitTo(*data,Save()) ;



  // P r i n t   f i t   r e s u l t s 
  // ---------------------------------

  // Summary printing: Basic info plus final values of floating fit parameters
  r->Print() ;

  // Verbose printing: Basic info, values of constant paramaters, initial and
  // final values of floating parameters, global correlations
  r->Print("v") ;



  // V i s u a l i z e   c o r r e l a t i o n   m a t r i x
  // -------------------------------------------------------

  // Construct 2D color plot of correlation matrix
  gStyle->SetOptStat(0) ;
  gStyle->SetPalette(1) ;
  TH2* hcorr = r->correlationHist() ;


  // Visualize ellipse corresponding to single correlation matrix element
  RooPlot* frame = new RooPlot(sigma1,sig1frac,0.45,0.60,0.65,0.90) ;
  frame->SetTitle("Covariance between sigma1 and sig1frac") ;
  r->plotOn(frame,sigma1,sig1frac,"ME12ABHV") ;



  // A c c e s s   f i t   r e s u l t   i n f o r m a t i o n
  // ---------------------------------------------------------

  // Access basic information
  cout << "EDM = " << r->edm() << endl ;
  cout << "-log(L) at minimum = " << r->minNll() << endl ;

  // Access list of final fit parameter values
  cout << "final value of floating parameters" << endl ;
  r->floatParsFinal().Print("s") ;

  // Access correlation matrix elements
  cout << "correlation between sig1frac and a0 is  " << r->correlation(sig1frac,a0) << endl ;
  cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac","mean") << endl ;

  // Extract covariance and correlation matrix as TMatrixDSym
  const TMatrixDSym& cor = r->correlationMatrix() ;
  const TMatrixDSym& cov = r->covarianceMatrix() ;

  // Print correlation, covariance matrix
  cout << "correlation matrix" << endl ;
  cor.Print() ;
  cout << "covariance matrix" << endl ;
  cov.Print() ;
  

  // P e r s i s t   f i t   r e s u l t   i n   r o o t   f i l e 
  // -------------------------------------------------------------
  
  // Open new ROOT file save save result 
  TFile f("rf607_fitresult.root","RECREATE") ;
  r->Write("rf607") ;
  f.Close() ;

  // In a clean ROOT session retrieve the persisted fit result as follows:
  // RooFitResult* r = gDirectory->Get("rf607") ;
 

  TCanvas* c = new TCanvas("rf607_fitresult","rf607_fitresult",800,400) ;
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; hcorr->GetYaxis()->SetTitleOffset(1.4) ; hcorr->Draw("colz") ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;

}
 rf607_fitresult.C:1
 rf607_fitresult.C:2
 rf607_fitresult.C:3
 rf607_fitresult.C:4
 rf607_fitresult.C:5
 rf607_fitresult.C:6
 rf607_fitresult.C:7
 rf607_fitresult.C:8
 rf607_fitresult.C:9
 rf607_fitresult.C:10
 rf607_fitresult.C:11
 rf607_fitresult.C:12
 rf607_fitresult.C:13
 rf607_fitresult.C:14
 rf607_fitresult.C:15
 rf607_fitresult.C:16
 rf607_fitresult.C:17
 rf607_fitresult.C:18
 rf607_fitresult.C:19
 rf607_fitresult.C:20
 rf607_fitresult.C:21
 rf607_fitresult.C:22
 rf607_fitresult.C:23
 rf607_fitresult.C:24
 rf607_fitresult.C:25
 rf607_fitresult.C:26
 rf607_fitresult.C:27
 rf607_fitresult.C:28
 rf607_fitresult.C:29
 rf607_fitresult.C:30
 rf607_fitresult.C:31
 rf607_fitresult.C:32
 rf607_fitresult.C:33
 rf607_fitresult.C:34
 rf607_fitresult.C:35
 rf607_fitresult.C:36
 rf607_fitresult.C:37
 rf607_fitresult.C:38
 rf607_fitresult.C:39
 rf607_fitresult.C:40
 rf607_fitresult.C:41
 rf607_fitresult.C:42
 rf607_fitresult.C:43
 rf607_fitresult.C:44
 rf607_fitresult.C:45
 rf607_fitresult.C:46
 rf607_fitresult.C:47
 rf607_fitresult.C:48
 rf607_fitresult.C:49
 rf607_fitresult.C:50
 rf607_fitresult.C:51
 rf607_fitresult.C:52
 rf607_fitresult.C:53
 rf607_fitresult.C:54
 rf607_fitresult.C:55
 rf607_fitresult.C:56
 rf607_fitresult.C:57
 rf607_fitresult.C:58
 rf607_fitresult.C:59
 rf607_fitresult.C:60
 rf607_fitresult.C:61
 rf607_fitresult.C:62
 rf607_fitresult.C:63
 rf607_fitresult.C:64
 rf607_fitresult.C:65
 rf607_fitresult.C:66
 rf607_fitresult.C:67
 rf607_fitresult.C:68
 rf607_fitresult.C:69
 rf607_fitresult.C:70
 rf607_fitresult.C:71
 rf607_fitresult.C:72
 rf607_fitresult.C:73
 rf607_fitresult.C:74
 rf607_fitresult.C:75
 rf607_fitresult.C:76
 rf607_fitresult.C:77
 rf607_fitresult.C:78
 rf607_fitresult.C:79
 rf607_fitresult.C:80
 rf607_fitresult.C:81
 rf607_fitresult.C:82
 rf607_fitresult.C:83
 rf607_fitresult.C:84
 rf607_fitresult.C:85
 rf607_fitresult.C:86
 rf607_fitresult.C:87
 rf607_fitresult.C:88
 rf607_fitresult.C:89
 rf607_fitresult.C:90
 rf607_fitresult.C:91
 rf607_fitresult.C:92
 rf607_fitresult.C:93
 rf607_fitresult.C:94
 rf607_fitresult.C:95
 rf607_fitresult.C:96
 rf607_fitresult.C:97
 rf607_fitresult.C:98
 rf607_fitresult.C:99
 rf607_fitresult.C:100
 rf607_fitresult.C:101
 rf607_fitresult.C:102
 rf607_fitresult.C:103
 rf607_fitresult.C:104
 rf607_fitresult.C:105
 rf607_fitresult.C:106
 rf607_fitresult.C:107
 rf607_fitresult.C:108
 rf607_fitresult.C:109
 rf607_fitresult.C:110
 rf607_fitresult.C:111
 rf607_fitresult.C:112
 rf607_fitresult.C:113
 rf607_fitresult.C:114
 rf607_fitresult.C:115
 rf607_fitresult.C:116
 rf607_fitresult.C:117
 rf607_fitresult.C:118
 rf607_fitresult.C:119
 rf607_fitresult.C:120
 rf607_fitresult.C:121
 rf607_fitresult.C:122
 rf607_fitresult.C:123
 rf607_fitresult.C:124
 rf607_fitresult.C:125
 rf607_fitresult.C:126
 rf607_fitresult.C:127
 rf607_fitresult.C:128
 rf607_fitresult.C:129
 rf607_fitresult.C:130
 rf607_fitresult.C:131
 rf607_fitresult.C:132
 rf607_fitresult.C:133
 rf607_fitresult.C:134
 rf607_fitresult.C:135
 rf607_fitresult.C:136
 rf607_fitresult.C:137
 rf607_fitresult.C:138
 rf607_fitresult.C:139
 rf607_fitresult.C:140
 rf607_fitresult.C:141
 rf607_fitresult.C:142
 rf607_fitresult.C:143
 rf607_fitresult.C:144
 rf607_fitresult.C:145
 rf607_fitresult.C:146
 rf607_fitresult.C:147
 rf607_fitresult.C:148