Logo ROOT  
Reference Guide
rf706_histpdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 ///
5 /// Special p.d.f.'s: histogram-based p.d.f.s and functions
6 ///
7 /// \macro_image
8 /// \macro_output
9 /// \macro_code
10 ///
11 /// \date 07/2008
12 /// \author Wouter Verkerke
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooPolynomial.h"
19 #include "RooHistPdf.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 using namespace RooFit;
24 
25 void rf706_histpdf()
26 {
27  // C r e a t e p d f f o r s a m p l i n g
28  // ---------------------------------------------
29 
30  RooRealVar x("x", "x", 0, 20);
31  RooPolynomial p("p", "p", x, RooArgList(RooConst(0.01), RooConst(-0.01), RooConst(0.0004)));
32 
33  // C r e a t e l o w s t a t s h i s t o g r a m
34  // ---------------------------------------------------
35 
36  // Sample 500 events from p
37  x.setBins(20);
38  RooDataSet *data1 = p.generate(x, 500);
39 
40  // Create a binned dataset with 20 bins and 500 events
41  RooDataHist *hist1 = data1->binnedClone();
42 
43  // Represent data in dh as pdf in x
44  RooHistPdf histpdf1("histpdf1", "histpdf1", x, *hist1, 0);
45 
46  // Plot unbinned data and histogram pdf overlaid
47  RooPlot *frame1 = x.frame(Title("Low statistics histogram pdf"), Bins(100));
48  data1->plotOn(frame1);
49  histpdf1.plotOn(frame1);
50 
51  // C r e a t e h i g h s t a t s h i s t o g r a m
52  // -----------------------------------------------------
53 
54  // Sample 100000 events from p
55  x.setBins(10);
56  RooDataSet *data2 = p.generate(x, 100000);
57 
58  // Create a binned dataset with 10 bins and 100K events
59  RooDataHist *hist2 = data2->binnedClone();
60 
61  // Represent data in dh as pdf in x, apply 2nd order interpolation
62  RooHistPdf histpdf2("histpdf2", "histpdf2", x, *hist2, 2);
63 
64  // Plot unbinned data and histogram pdf overlaid
65  RooPlot *frame2 = x.frame(Title("High stats histogram pdf with interpolation"), Bins(100));
66  data2->plotOn(frame2);
67  histpdf2.plotOn(frame2);
68 
69  TCanvas *c = new TCanvas("rf706_histpdf", "rf706_histpdf", 800, 400);
70  c->Divide(2);
71  c->cd(1);
72  gPad->SetLeftMargin(0.15);
73  frame1->GetYaxis()->SetTitleOffset(1.4);
74  frame1->Draw();
75  c->cd(2);
76  gPad->SetLeftMargin(0.15);
77  frame2->GetYaxis()->SetTitleOffset(1.8);
78  frame2->Draw();
79 }
c
#define c(i)
Definition: RSha256.hxx:119
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
RooFit::Bins
RooCmdArg Bins(Int_t nbin)
Definition: RooGlobalFunc.cxx:174
rf706_histpdf
Definition: rf706_histpdf.py:1
RooArgList
Definition: RooArgList.h:21
RooGaussian.h
x
Double_t x[n]
Definition: legend1.C:17
TCanvas.h
RooDataSet.h
RooPolynomial.h
RooDataHist
Definition: RooDataHist.h:39
RooFit
Definition: RooCFunction1Binding.h:29
RooPolynomial
Definition: RooPolynomial.h:28
RooAbsData::plotOn
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Definition: RooAbsData.cxx:547
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1256
RooPlot
Definition: RooPlot.h:44
RooRealVar.h
RooHistPdf
Definition: RooHistPdf.h:29
RooConstVar.h
RooHistPdf.h
TCanvas
Definition: TCanvas.h:23
TAxis.h
RooDataSet::binnedClone
RooDataHist * binnedClone(const char *newName=0, const char *newTitle=0) const
Return binned clone of this dataset.
Definition: RooDataSet.cxx:972
gPad
#define gPad
Definition: TVirtualPad.h:287
RooDataSet
Definition: RooDataSet.h:33
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
Definition: RooRealVar.h:35
RooFit::Title
RooCmdArg Title(const char *name)
Definition: RooGlobalFunc.cxx:173
RooFit::RooConst
RooConstVar & RooConst(Double_t val)
Definition: RooGlobalFunc.cxx:341