ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605
// 
// Working with the profile likelihood estimator
//
//
//
// 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 "RooMinuit.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf605_profilell()
{
  // 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,-10,10) ;
  RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
  RooGaussian g1("g1","g1",x,mean,sigma_g1) ;

  RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.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) ;
  


  // C o n s t r u c t   p l a i n   l i k e l i h o o d
  // ---------------------------------------------------

  // Construct unbinned likelihood
  RooAbsReal* nll = model.createNLL(*data,NumCPU(2)) ;

  // Minimize likelihood w.r.t all parameters before making plots
  RooMinuit(*nll).migrad() ;

  // Plot likelihood scan frac 
  RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
  nll->plotOn(frame1,ShiftToZero()) ;

  // Plot likelihood scan in sigma_g2
  RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
  nll->plotOn(frame2,ShiftToZero()) ;



  // C o n s t r u c t   p r o f i l e   l i k e l i h o o d   i n   f r a c
  // -----------------------------------------------------------------------

  // The profile likelihood estimator on nll for frac will minimize nll w.r.t
  // all floating parameters except frac for each evaluation

  RooAbsReal* pll_frac = nll->createProfile(frac) ;

  // Plot the profile likelihood in frac
  pll_frac->plotOn(frame1,LineColor(kRed)) ;

  // Adjust frame maximum for visual clarity
  frame1->SetMinimum(0) ;
  frame1->SetMaximum(3) ;



  // C o n s t r u c t   p r o f i l e   l i k e l i h o o d   i n   s i g m a _ g 2 
  // -------------------------------------------------------------------------------

  // The profile likelihood estimator on nll for sigma_g2 will minimize nll
  // w.r.t all floating parameters except sigma_g2 for each evaluation
  RooAbsReal* pll_sigmag2 = nll->createProfile(sigma_g2) ;

  // Plot the profile likelihood in sigma_g2
  pll_sigmag2->plotOn(frame2,LineColor(kRed)) ;

  // Adjust frame maximum for visual clarity
  frame2->SetMinimum(0) ;
  frame2->SetMaximum(3) ;



  // Make canvas and draw RooPlots
  TCanvas *c = new TCanvas("rf605_profilell","rf605_profilell",800, 400);
  c->Divide(2);
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;

  delete pll_frac ;
  delete pll_sigmag2 ;
  delete nll ;
}

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