ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #701
// 
// Unbinned maximum likelihood fit of an efficiency eff(x) function to 
// a dataset D(x,cut), where cut is a category encoding a selection, of which
// the efficiency as function of x should be described by eff(x)
//
// 07/2008 - Wouter Verkerke 
//
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooFormulaVar.h"
#include "RooProdPdf.h"
#include "RooEfficiency.h"
#include "RooPolynomial.h"
#include "RooCategory.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf701_efficiencyfit()
{
  // C o n s t r u c t   e f f i c i e n c y   f u n c t i o n   e ( x ) 
  // -------------------------------------------------------------------

  // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
  RooRealVar x("x","x",-10,10) ;

  // Efficiency function eff(x;a,b) 
  RooRealVar a("a","a",0.4,0,1) ;
  RooRealVar b("b","b",5) ;
  RooRealVar c("c","c",-1,-10,10) ;
  RooFormulaVar effFunc("effFunc","(1-a)+a*cos((x-c)/b)",RooArgList(a,b,c,x)) ;



  // C o n s t r u c t   c o n d i t i o n a l    e f f i c i e n c y   p d f   E ( c u t | x ) 
  // ------------------------------------------------------------------------------------------

  // Acceptance state cut (1 or 0)
  RooCategory cut("cut","cutr") ;
  cut.defineType("accept",1) ;
  cut.defineType("reject",0) ;

  // Construct efficiency p.d.f eff(cut|x)
  RooEfficiency effPdf("effPdf","effPdf",effFunc,cut,"accept") ;



  // G e n e r a t e   d a t a   ( x ,   c u t )   f r o m   a   t o y   m o d e l 
  // -----------------------------------------------------------------------------

  // Construct global shape p.d.f shape(x) and product model(x,cut) = eff(cut|x)*shape(x) 
  // (These are _only_ needed to generate some toy MC here to be used later)
  RooPolynomial shapePdf("shapePdf","shapePdf",x,RooConst(-0.095)) ;
  RooProdPdf model("model","model",shapePdf,Conditional(effPdf,cut)) ;

  // Generate some toy data from model
  RooDataSet* data = model.generate(RooArgSet(x,cut),10000) ;



  // F i t   c o n d i t i o n a l   e f f i c i e n c y   p d f   t o   d a t a 
  // --------------------------------------------------------------------------

  // Fit conditional efficiency p.d.f to data
  effPdf.fitTo(*data,ConditionalObservables(x)) ;



  // P l o t   f i t t e d ,   d a t a   e f f i c i e n c y  
  // --------------------------------------------------------

  // Plot distribution of all events and accepted fraction of events on frame
  RooPlot* frame1 = x.frame(Bins(20),Title("Data (all, accepted)")) ;
  data->plotOn(frame1) ;
  data->plotOn(frame1,Cut("cut==cut::accept"),MarkerColor(kRed),LineColor(kRed)) ;

  // Plot accept/reject efficiency on data overlay fitted efficiency curve
  RooPlot* frame2 = x.frame(Bins(20),Title("Fitted efficiency")) ;
  data->plotOn(frame2,Efficiency(cut)) ; // needs ROOT version >= 5.21
  effFunc.plotOn(frame2,LineColor(kRed)) ;



  // Draw all frames on a canvas
  TCanvas* ca = new TCanvas("rf701_efficiency","rf701_efficiency",800,400) ;
  ca->Divide(2) ;
  ca->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.6) ; frame1->Draw() ;
  ca->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
  
 
}
 rf701_efficiencyfit.C:1
 rf701_efficiencyfit.C:2
 rf701_efficiencyfit.C:3
 rf701_efficiencyfit.C:4
 rf701_efficiencyfit.C:5
 rf701_efficiencyfit.C:6
 rf701_efficiencyfit.C:7
 rf701_efficiencyfit.C:8
 rf701_efficiencyfit.C:9
 rf701_efficiencyfit.C:10
 rf701_efficiencyfit.C:11
 rf701_efficiencyfit.C:12
 rf701_efficiencyfit.C:13
 rf701_efficiencyfit.C:14
 rf701_efficiencyfit.C:15
 rf701_efficiencyfit.C:16
 rf701_efficiencyfit.C:17
 rf701_efficiencyfit.C:18
 rf701_efficiencyfit.C:19
 rf701_efficiencyfit.C:20
 rf701_efficiencyfit.C:21
 rf701_efficiencyfit.C:22
 rf701_efficiencyfit.C:23
 rf701_efficiencyfit.C:24
 rf701_efficiencyfit.C:25
 rf701_efficiencyfit.C:26
 rf701_efficiencyfit.C:27
 rf701_efficiencyfit.C:28
 rf701_efficiencyfit.C:29
 rf701_efficiencyfit.C:30
 rf701_efficiencyfit.C:31
 rf701_efficiencyfit.C:32
 rf701_efficiencyfit.C:33
 rf701_efficiencyfit.C:34
 rf701_efficiencyfit.C:35
 rf701_efficiencyfit.C:36
 rf701_efficiencyfit.C:37
 rf701_efficiencyfit.C:38
 rf701_efficiencyfit.C:39
 rf701_efficiencyfit.C:40
 rf701_efficiencyfit.C:41
 rf701_efficiencyfit.C:42
 rf701_efficiencyfit.C:43
 rf701_efficiencyfit.C:44
 rf701_efficiencyfit.C:45
 rf701_efficiencyfit.C:46
 rf701_efficiencyfit.C:47
 rf701_efficiencyfit.C:48
 rf701_efficiencyfit.C:49
 rf701_efficiencyfit.C:50
 rf701_efficiencyfit.C:51
 rf701_efficiencyfit.C:52
 rf701_efficiencyfit.C:53
 rf701_efficiencyfit.C:54
 rf701_efficiencyfit.C:55
 rf701_efficiencyfit.C:56
 rf701_efficiencyfit.C:57
 rf701_efficiencyfit.C:58
 rf701_efficiencyfit.C:59
 rf701_efficiencyfit.C:60
 rf701_efficiencyfit.C:61
 rf701_efficiencyfit.C:62
 rf701_efficiencyfit.C:63
 rf701_efficiencyfit.C:64
 rf701_efficiencyfit.C:65
 rf701_efficiencyfit.C:66
 rf701_efficiencyfit.C:67
 rf701_efficiencyfit.C:68
 rf701_efficiencyfit.C:69
 rf701_efficiencyfit.C:70
 rf701_efficiencyfit.C:71
 rf701_efficiencyfit.C:72
 rf701_efficiencyfit.C:73
 rf701_efficiencyfit.C:74
 rf701_efficiencyfit.C:75
 rf701_efficiencyfit.C:76
 rf701_efficiencyfit.C:77
 rf701_efficiencyfit.C:78
 rf701_efficiencyfit.C:79
 rf701_efficiencyfit.C:80
 rf701_efficiencyfit.C:81
 rf701_efficiencyfit.C:82
 rf701_efficiencyfit.C:83
 rf701_efficiencyfit.C:84
 rf701_efficiencyfit.C:85
 rf701_efficiencyfit.C:86
 rf701_efficiencyfit.C:87
 rf701_efficiencyfit.C:88
 rf701_efficiencyfit.C:89
 rf701_efficiencyfit.C:90
 rf701_efficiencyfit.C:91
 rf701_efficiencyfit.C:92
 rf701_efficiencyfit.C:93
 rf701_efficiencyfit.C:94
 rf701_efficiencyfit.C:95
 rf701_efficiencyfit.C:96
 rf701_efficiencyfit.C:97
 rf701_efficiencyfit.C:98
 rf701_efficiencyfit.C:99
 rf701_efficiencyfit.C:100
 rf701_efficiencyfit.C:101
 rf701_efficiencyfit.C:102
 rf701_efficiencyfit.C:103
 rf701_efficiencyfit.C:104