Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf706_histpdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Special p.d.f.'s: histogram-based p.d.f.s and functions
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooPolynomial.h"
16 #include "RooHistPdf.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "RooPlot.h"
20 using namespace RooFit;
21 
22 void rf706_histpdf()
23 {
24  // C r e a t e p d f f o r s a m p l i n g
25  // ---------------------------------------------
26 
27  RooRealVar x("x", "x", 0, 20);
28  RooPolynomial p("p", "p", x, RooArgList(RooConst(0.01), RooConst(-0.01), RooConst(0.0004)));
29 
30  // C r e a t e l o w s t a t s h i s t o g r a m
31  // ---------------------------------------------------
32 
33  // Sample 500 events from p
34  x.setBins(20);
35  RooDataSet *data1 = p.generate(x, 500);
36 
37  // Create a binned dataset with 20 bins and 500 events
38  RooDataHist *hist1 = data1->binnedClone();
39 
40  // Represent data in dh as pdf in x
41  RooHistPdf histpdf1("histpdf1", "histpdf1", x, *hist1, 0);
42 
43  // Plot unbinned data and histogram pdf overlaid
44  RooPlot *frame1 = x.frame(Title("Low statistics histogram pdf"), Bins(100));
45  data1->plotOn(frame1);
46  histpdf1.plotOn(frame1);
47 
48  // C r e a t e h i g h s t a t s h i s t o g r a m
49  // -----------------------------------------------------
50 
51  // Sample 100000 events from p
52  x.setBins(10);
53  RooDataSet *data2 = p.generate(x, 100000);
54 
55  // Create a binned dataset with 10 bins and 100K events
56  RooDataHist *hist2 = data2->binnedClone();
57 
58  // Represent data in dh as pdf in x, apply 2nd order interpolation
59  RooHistPdf histpdf2("histpdf2", "histpdf2", x, *hist2, 2);
60 
61  // Plot unbinned data and histogram pdf overlaid
62  RooPlot *frame2 = x.frame(Title("High stats histogram pdf with interpolation"), Bins(100));
63  data2->plotOn(frame2);
64  histpdf2.plotOn(frame2);
65 
66  TCanvas *c = new TCanvas("rf706_histpdf", "rf706_histpdf", 800, 400);
67  c->Divide(2);
68  c->cd(1);
69  gPad->SetLeftMargin(0.15);
70  frame1->GetYaxis()->SetTitleOffset(1.4);
71  frame1->Draw();
72  c->cd(2);
73  gPad->SetLeftMargin(0.15);
74  frame2->GetYaxis()->SetTitleOffset(1.8);
75  frame2->Draw();
76 }