rf706_histpdf

Special pdf's: histogram-based pdfs and functions

Author: Wouter Verkerke
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, November 30, 2022 at 11:24 AM.

In [1]:
%%cpp -d
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooPolynomial.h"
#include "RooHistPdf.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;

Create pdf for sampling

In [2]:
RooRealVar x("x", "x", 0, 20);
RooPolynomial p("p", "p", x, RooArgList(0.01, -0.01, 0.0004));

Create low stats histogram

Sample 500 events from p

In [3]:
x.setBins(20);
RooDataSet *data1 = p.generate(x, 500);

Create a binned dataset with 20 bins and 500 events

In [4]:
RooDataHist *hist1 = data1->binnedClone();

Represent data in dh as pdf in x

In [5]:
RooHistPdf histpdf1("histpdf1", "histpdf1", x, *hist1, 0);

Plot unbinned data and histogram pdf overlaid

In [6]:
RooPlot *frame1 = x.frame(Title("Low statistics histogram pdf"), Bins(100));
data1->plotOn(frame1);
histpdf1.plotOn(frame1);

Create high stats histogram

Sample 100000 events from p

In [7]:
x.setBins(10);
RooDataSet *data2 = p.generate(x, 100000);

Create a binned dataset with 10 bins and 100K events

In [8]:
RooDataHist *hist2 = data2->binnedClone();

Represent data in dh as pdf in x, apply 2nd order interpolation

In [9]:
RooHistPdf histpdf2("histpdf2", "histpdf2", x, *hist2, 2);

Plot unbinned data and histogram pdf overlaid

In [10]:
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();

Draw all canvases

In [11]:
%jsroot on
gROOT->GetListOfCanvases()->Draw()