ROOT logo

From $ROOTSYS/tutorials/roofit/rf702_efficiencyfit_2D.C

//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #702
// 
// Unbinned maximum likelihood fit of an efficiency eff(x) function to 
// a dataset D(x,cut), where cut is a category encoding a selection whose 
// efficiency as function of x should be described by eff(x)
//
//
/////////////////////////////////////////////////////////////////////////

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


void rf702_efficiencyfit_2D(Bool_t flat=kFALSE)
{
  // 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 , y ) 
  // -----------------------------------------------------------------------

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

  // Efficiency function eff(x;a,b) 
  RooRealVar ax("ax","ay",0.6,0,1) ;
  RooRealVar bx("bx","by",5) ;
  RooRealVar cx("cx","cy",-1,-10,10) ;

  RooRealVar ay("ay","ay",0.2,0,1) ;
  RooRealVar by("by","by",5) ;
  RooRealVar cy("cy","cy",-1,-10,10) ;

  RooFormulaVar effFunc("effFunc","((1-ax)+ax*cos((x-cx)/bx))*((1-ay)+ay*cos((y-cy)/by))",RooArgList(ax,bx,cx,x,ay,by,cy,y)) ; 

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



  // 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 , y ) 
  // ---------------------------------------------------------------------------------------------

  // 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 , y , 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 shapePdfX("shapePdfX","shapePdfX",x,RooConst(flat?0:-0.095)) ;
  RooPolynomial shapePdfY("shapePdfY","shapePdfY",y,RooConst(flat?0:+0.095)) ;
  RooProdPdf shapePdf("shapePdf","shapePdf",RooArgSet(shapePdfX,shapePdfY)) ;
  RooProdPdf model("model","model",shapePdf,Conditional(effPdf,cut)) ;

  // Generate some toy data from model
  RooDataSet* data = model.generate(RooArgSet(x,y,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(RooArgSet(x,y))) ;



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

  // Make 2D histograms of all data, selected data and efficiency function
  TH1* hh_data_all = data->createHistogram("hh_data_all",x,Binning(8),YVar(y,Binning(8))) ;
  TH1* hh_data_sel = data->createHistogram("hh_data_sel",x,Binning(8),YVar(y,Binning(8)),Cut("cut==cut::accept")) ;
  TH1* hh_eff      = effFunc.createHistogram("hh_eff",x,Binning(50),YVar(y,Binning(50))) ;

  // Some adjustsment for good visualization
  hh_data_all->SetMinimum(0) ;
  hh_data_sel->SetMinimum(0) ;
  hh_eff->SetMinimum(0) ;
  hh_eff->SetLineColor(kBlue) ;



  // Draw all frames on a canvas
  TCanvas* ca = new TCanvas("rf702_efficiency_2D","rf702_efficiency_2D",1200,400) ;
  ca->Divide(3) ;
  ca->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_data_all->GetZaxis()->SetTitleOffset(1.8) ; hh_data_all->Draw("lego") ;
  ca->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_data_sel->GetZaxis()->SetTitleOffset(1.8) ; hh_data_sel->Draw("lego") ;
  ca->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_eff->GetZaxis()->SetTitleOffset(1.8) ; hh_eff->Draw("surf") ;

  return ;

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