Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf706_histpdf.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Special pdf's: histogram-based pdfs and functions
5///
6/// \macro_image
7/// \macro_code
8/// \macro_output
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooPolynomial.h"
17#include "RooHistPdf.h"
18#include "TCanvas.h"
19#include "TAxis.h"
20#include "RooPlot.h"
21using namespace RooFit;
22
23void rf706_histpdf()
24{
25 // C r e a t e p d f f o r s a m p l i n g
26 // ---------------------------------------------
27
28 RooRealVar x("x", "x", 0, 20);
29 RooPolynomial p("p", "p", x, RooArgList(0.01, -0.01, 0.0004));
30
31 // C r e a t e l o w s t a t s h i s t o g r a m
32 // ---------------------------------------------------
33
34 // Sample 500 events from p
35 x.setBins(20);
36 std::unique_ptr<RooDataSet> data1{p.generate(x, 500)};
37
38 // Create a binned dataset with 20 bins and 500 events
39 std::unique_ptr<RooDataHist> hist1{data1->binnedClone()};
40
41 // Represent data in dh as pdf in x
42 RooHistPdf histpdf1("histpdf1", "histpdf1", x, *hist1, 0);
43
44 // Plot unbinned data and histogram pdf overlaid
45 RooPlot *frame1 = x.frame(Title("Low statistics histogram pdf"), Bins(100));
46 data1->plotOn(frame1);
47 histpdf1.plotOn(frame1);
48
49 // C r e a t e h i g h s t a t s h i s t o g r a m
50 // -----------------------------------------------------
51
52 // Sample 100000 events from p
53 x.setBins(10);
54 std::unique_ptr<RooDataSet> data2{p.generate(x, 100000)};
55
56 // Create a binned dataset with 10 bins and 100K events
57 std::unique_ptr<RooDataHist> hist2{data2->binnedClone()};
58
59 // Represent data in dh as pdf in x, apply 2nd order interpolation
60 RooHistPdf histpdf2("histpdf2", "histpdf2", x, *hist2, 2);
61
62 // Plot unbinned data and histogram pdf overlaid
63 RooPlot *frame2 = x.frame(Title("High stats histogram pdf with interpolation"), Bins(100));
64 data2->plotOn(frame2);
65 histpdf2.plotOn(frame2);
66
67 TCanvas *c = new TCanvas("rf706_histpdf", "rf706_histpdf", 800, 400);
68 c->Divide(2);
69 c->cd(1);
70 gPad->SetLeftMargin(0.15);
71 frame1->GetYaxis()->SetTitleOffset(1.4);
72 frame1->Draw();
73 c->cd(2);
74 gPad->SetLeftMargin(0.15);
75 frame2->GetYaxis()->SetTitleOffset(1.8);
76 frame2->Draw();
77}
#define c(i)
Definition RSha256.hxx:101
winID h TVirtualViewer3D TVirtualGLPainter p
#define gPad
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
A propability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
RooPolynomial implements a polynomial p.d.f of the form.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68