ROOT logo

From $ROOTSYS/tutorials/roofit/rf608_fitresultaspdf.C

//////////////////////////////////////////////////////////////////////////
//
// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #608
// 
// Representing the parabolic approximation of the fit as
// a multi-variate Gaussian on the parameters of the fitted p.d.f.
//
//
// 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 "TH3.h"

using namespace RooFit ;


void rf608_fitresultaspdf()
{
  // C r e a t e   m o d e l   a n d   d a t a s e t 
  // -----------------------------------------------

  // Observable
  RooRealVar x("x","x",-20,20) ;

  // Model (intentional strong correlations)
  RooRealVar mean("mean","mean of g1 and g2",0,-1,1) ;
  RooRealVar sigma_g1("sigma_g1","width of g1",2) ; 
  RooGaussian g1("g1","g1",x,mean,sigma_g1) ;

  RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,5.0) ;
  RooGaussian g2("g2","g2",x,mean,sigma_g2) ;

  RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
  RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;

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


  // F i t   m o d e l   t o   d a t a 
  // ----------------------------------

  RooFitResult* r = model.fitTo(*data,Save()) ;


  // C r e a t e M V   G a u s s i a n   p d f   o f   f i t t e d    p a r a m e t e r s
  // ------------------------------------------------------------------------------------

  RooAbsPdf* parabPdf = r->createHessePdf(RooArgSet(frac,mean,sigma_g2)) ;


  // S o m e   e x e c e r c i s e s   w i t h   t h e   p a r a m e t e r   p d f 
  // -----------------------------------------------------------------------------

  // Generate 100K points in the parameter space, sampled from the MVGaussian p.d.f.
  RooDataSet* d = parabPdf->generate(RooArgSet(mean,sigma_g2,frac),100000) ;


  // Sample a 3-D histogram of the p.d.f. to be visualized as an error ellipsoid using the GLISO draw option
  TH3* hh_3d = (TH3*) parabPdf->createHistogram("mean,sigma_g2,frac",25,25,25) ;
  hh_3d->SetFillColor(kBlue) ;  


  // Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s 
  // The integrations corresponding to these projections are performed analytically
  // by the MV Gaussian p.d.f.
  RooAbsPdf* pdf_sigmag2_frac = parabPdf->createProjection(mean) ;
  RooAbsPdf* pdf_mean_frac    = parabPdf->createProjection(sigma_g2) ;
  RooAbsPdf* pdf_mean_sigmag2 = parabPdf->createProjection(frac) ;


  // Make 2D plots of the 3 two-dimensional p.d.f. projections
  TH2* hh_sigmag2_frac = (TH2*) pdf_sigmag2_frac->createHistogram("sigma_g2,frac",50,50) ;
  TH2* hh_mean_frac    = (TH2*) pdf_mean_frac->createHistogram("mean,frac",50,50) ;
  TH2* hh_mean_sigmag2 = (TH2*) pdf_mean_sigmag2->createHistogram("mean,sigma_g2",50,50) ;
  hh_mean_frac->SetLineColor(kBlue) ;
  hh_sigmag2_frac->SetLineColor(kBlue) ;
  hh_mean_sigmag2->SetLineColor(kBlue) ;


  // Draw the 'sigar'
  gStyle->SetCanvasPreferGL(true);
  gStyle->SetPalette(1) ;
  new TCanvas("rf608_fitresultaspdf_1","rf608_fitresultaspdf_1",600,600) ;
  hh_3d->Draw("gliso") ; 

  // Draw the 2D projections of the 3D p.d.f.
  TCanvas* c2 = new TCanvas("rf608_fitresultaspdf_2","rf608_fitresultaspdf_2",900,600) ;
  c2->Divide(3,2) ;
  c2->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_sigmag2->Draw("surf3") ; 
  c2->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_sigmag2_frac->Draw("surf3") ; 
  c2->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_mean_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_frac->Draw("surf3") ; 

  // Draw the distributions of parameter points sampled from the p.d.f.
  TH1* tmp1 = d->createHistogram("mean,sigma_g2",50,50) ;
  TH1* tmp2 = d->createHistogram("sigma_g2,frac",50,50) ;
  TH1* tmp3 = d->createHistogram("mean,frac",50,50) ;

  c2->cd(4) ; gPad->SetLeftMargin(0.15) ; tmp1->GetZaxis()->SetTitleOffset(1.4) ; tmp1->Draw("lego3") ;
  c2->cd(5) ; gPad->SetLeftMargin(0.15) ; tmp2->GetZaxis()->SetTitleOffset(1.4) ; tmp2->Draw("lego3") ;
  c2->cd(6) ; gPad->SetLeftMargin(0.15) ; tmp3->GetZaxis()->SetTitleOffset(1.4) ; tmp3->Draw("lego3") ;

}

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