ROOT logo

From $ROOTSYS/tutorials/roofit/rf703_effpdfprod.C

//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #703
// 
// Using a product of an (acceptance) efficiency and a p.d.f as 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 "RooExponential.h"
#include "RooEffProd.h"
#include "RooFormulaVar.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf703_effpdfprod()
{
  // D e f i n e   o b s e r v a b l e s   a n d   d e c a y   p d f
  // ---------------------------------------------------------------

  // Declare observables
  RooRealVar t("t","t",0,5) ;

  // Make pdf
  RooRealVar tau("tau","tau",-1.54,-4,-0.1) ;
  RooExponential model("model","model",t,tau) ;



  // D e f i n e   e f f i c i e n c y   f u n c t i o n
  // ---------------------------------------------------
  
  // Use error function to simulate turn-on slope
  RooFormulaVar eff("eff","0.5*(TMath::Erf((t-1)/0.5)+1)",t) ;



  // D e f i n e   d e c a y   p d f   w i t h   e f f i c i e n c y 
  // ---------------------------------------------------------------

  // Multiply pdf(t) with efficiency in t
  RooEffProd modelEff("modelEff","model with efficiency",model,eff) ;

  

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

  RooPlot* frame1 = t.frame(Title("Efficiency")) ;
  eff.plotOn(frame1,LineColor(kRed)) ;

  RooPlot* frame2 = t.frame(Title("Pdf with and without efficiency")) ;

  model.plotOn(frame2,LineStyle(kDashed)) ;
  modelEff.plotOn(frame2) ;



  // G e n e r a t e   t o y   d a t a ,   f i t   m o d e l E f f   t o   d a t a
  // ------------------------------------------------------------------------------

  // Generate events. If the input pdf has an internal generator, the internal generator
  // is used and an accept/reject sampling on the efficiency is applied. 
  RooDataSet* data = modelEff.generate(t,10000) ;

  // Fit pdf. The normalization integral is calculated numerically. 
  modelEff.fitTo(*data) ;

  // Plot generated data and overlay fitted pdf
  RooPlot* frame3 = t.frame(Title("Fitted pdf with efficiency")) ;
  data->plotOn(frame3) ;
  modelEff.plotOn(frame3) ;

  
  TCanvas* c = new TCanvas("rf703_effpdfprod","rf703_effpdfprod",1200,400) ;
  c->Divide(3) ;
  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.6) ; frame2->Draw() ;
  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;

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