Logo ROOT   6.18/05
Reference Guide
rf105_funcbinding.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Basic functionality: binding ROOT math functions as RooFit functions and pdfs
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 "TCanvas.h"
15#include "TAxis.h"
16#include "RooPlot.h"
17#include "TMath.h"
18#include "TF1.h"
19#include "Math/DistFunc.h"
20#include "RooTFnBinding.h"
21
22using namespace RooFit;
23
24void rf105_funcbinding()
25{
26
27 // B i n d T M a t h : : E r f C f u n c t i o n
28 // ---------------------------------------------------
29
30 // Bind one-dimensional TMath::Erf function as RooAbsReal function
31 RooRealVar x("x", "x", -3, 3);
32 RooAbsReal *errorFunc = bindFunction("erf", TMath::Erf, x);
33
34 // Print erf definition
35 errorFunc->Print();
36
37 // Plot erf on frame
38 RooPlot *frame1 = x.frame(Title("TMath::Erf bound as RooFit function"));
39 errorFunc->plotOn(frame1);
40
41 // B i n d R O O T : : M a t h : : b e t a _ p d f C f u n c t i o n
42 // -----------------------------------------------------------------------
43
44 // Bind pdf ROOT::Math::Beta with three variables as RooAbsPdf function
45 RooRealVar x2("x2", "x2", 0, 0.999);
46 RooRealVar a("a", "a", 5, 0, 10);
47 RooRealVar b("b", "b", 2, 0, 10);
49
50 // Perf beta definition
51 beta->Print();
52
53 // Generate some events and fit
54 RooDataSet *data = beta->generate(x2, 10000);
55 beta->fitTo(*data);
56
57 // Plot data and pdf on frame
58 RooPlot *frame2 = x2.frame(Title("ROOT::Math::Beta bound as RooFit pdf"));
59 data->plotOn(frame2);
60 beta->plotOn(frame2);
61
62 // B i n d R O O T T F 1 a s R o o F i t f u n c t i o n
63 // ---------------------------------------------------------------
64
65 // Create a ROOT TF1 function
66 TF1 *fa1 = new TF1("fa1", "sin(x)/x", 0, 10);
67
68 // Create an observable
69 RooRealVar x3("x3", "x3", 0.01, 20);
70
71 // Create binding of TF1 object to above observable
72 RooAbsReal *rfa1 = bindFunction(fa1, x3);
73
74 // Print rfa1 definition
75 rfa1->Print();
76
77 // Make plot frame in observable, plot TF1 binding function
78 RooPlot *frame3 = x3.frame(Title("TF1 bound as RooFit function"));
79 rfa1->plotOn(frame3);
80
81 TCanvas *c = new TCanvas("rf105_funcbinding", "rf105_funcbinding", 1200, 400);
82 c->Divide(3);
83 c->cd(1);
84 gPad->SetLeftMargin(0.15);
85 frame1->GetYaxis()->SetTitleOffset(1.6);
86 frame1->Draw();
87 c->cd(2);
88 gPad->SetLeftMargin(0.15);
89 frame2->GetYaxis()->SetTitleOffset(1.6);
90 frame2->Draw();
91 c->cd(3);
92 gPad->SetLeftMargin(0.15);
93 frame3->GetYaxis()->SetTitleOffset(1.6);
94 frame3->Draw();
95}
#define b(i)
Definition: RSha256.hxx:100
#define c(i)
Definition: RSha256.hxx:101
static const double x2[5]
static const double x3[11]
#define gPad
Definition: TVirtualPad.h:286
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:272
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
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:552
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
1-Dim function class
Definition: TF1.h:211
double beta_pdf(double x, double a, double b)
Probability density function of the beta distribution.
double beta(double x, double y)
Calculates the beta function.
Double_t x[n]
Definition: legend1.C:17
Template specialisation used in RooAbsArg:
RooAbsPdf * bindPdf(const char *name, CFUNCD1D func, RooAbsReal &x)
RooAbsReal * bindFunction(const char *name, CFUNCD1D func, RooAbsReal &x)
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition: TMath.cxx:184
const char * Title
Definition: TXMLSetup.cxx:67
auto * a
Definition: textangle.C:12