ROOT logo

From $ROOTSYS/tutorials/roofit/rf706_histpdf.C

//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #706
// 
// Histogram based p.d.f.s and functions
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

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


void rf706_histpdf()
{
  // C r e a t e   p d f   f o r   s a m p l i n g 
  // ---------------------------------------------

  RooRealVar x("x","x",0,20) ;
  RooPolynomial p("p","p",x,RooArgList(RooConst(0.01),RooConst(-0.01),RooConst(0.0004))) ;



  // C r e a t e   l o w   s t a t s   h i s t o g r a m
  // ---------------------------------------------------

  // Sample 500 events from p
  x.setBins(20) ;
  RooDataSet* data1 = p.generate(x,500) ;
  
  // Create a binned dataset with 20 bins and 500 events
  RooDataHist* hist1 = data1->binnedClone() ;

  // Represent data in dh as pdf in x
  RooHistPdf histpdf1("histpdf1","histpdf1",x,*hist1,0) ;

  // Plot unbinned data and histogram pdf overlaid
  RooPlot* frame1 = x.frame(Title("Low statistics histogram pdf"),Bins(100)) ;
  data1->plotOn(frame1) ;
  histpdf1.plotOn(frame1) ;  
  

  // C r e a t e   h i g h   s t a t s   h i s t o g r a m
  // -----------------------------------------------------

  // Sample 100000 events from p
  x.setBins(10) ;
  RooDataSet* data2 = p.generate(x,100000) ;

  // Create a binned dataset with 10 bins and 100K events
  RooDataHist* hist2 = data2->binnedClone() ;

  // Represent data in dh as pdf in x, apply 2nd order interpolation  
  RooHistPdf histpdf2("histpdf2","histpdf2",x,*hist2,2) ;

  // Plot unbinned data and histogram pdf overlaid
  RooPlot* frame2 = x.frame(Title("High stats histogram pdf with interpolation"),Bins(100)) ;
  data2->plotOn(frame2) ;
  histpdf2.plotOn(frame2) ;  


  TCanvas* c = new TCanvas("rf706_histpdf","rf706_histpdf",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.8) ; frame2->Draw() ;


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