Logo ROOT   6.10/09
Reference Guide
rf105_funcbinding.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'BASIC FUNCTIONALITY' RooFit tutorial macro #105
5 ///
6 /// Demonstration of binding ROOT Math functions as RooFit functions
7 /// and pdfs
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "RooPlot.h"
20 #include "TMath.h"
21 #include "TF1.h"
22 #include "Math/DistFunc.h"
23 #include "RooTFnBinding.h"
24 
25 using namespace RooFit;
26 
27 
28 void rf105_funcbinding()
29 {
30 
31  // B i n d T M a t h : : E r f C f u n c t i o n
32  // ---------------------------------------------------
33 
34  // Bind one-dimensional TMath::Erf function as RooAbsReal function
35  RooRealVar x("x","x",-3,3) ;
36  RooAbsReal* errorFunc = bindFunction("erf",TMath::Erf,x) ;
37 
38  // Print erf definition
39  errorFunc->Print() ;
40 
41  // Plot erf on frame
42  RooPlot* frame1 = x.frame(Title("TMath::Erf bound as RooFit function")) ;
43  errorFunc->plotOn(frame1) ;
44 
45 
46  // 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
47  // -----------------------------------------------------------------------
48 
49  // Bind pdf ROOT::Math::Beta with three variables as RooAbsPdf function
50  RooRealVar x2("x2","x2",0,0.999) ;
51  RooRealVar a("a","a",5,0,10) ;
52  RooRealVar b("b","b",2,0,10) ;
54 
55  // Perf beta definition
56  beta->Print() ;
57 
58  // Generate some events and fit
59  RooDataSet* data = beta->generate(x2,10000) ;
60  beta->fitTo(*data) ;
61 
62  // Plot data and pdf on frame
63  RooPlot* frame2 = x2.frame(Title("ROOT::Math::Beta bound as RooFit pdf")) ;
64  data->plotOn(frame2) ;
65  beta->plotOn(frame2) ;
66 
67 
68 
69  // 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
70  // ---------------------------------------------------------------
71 
72  // Create a ROOT TF1 function
73  TF1 *fa1 = new TF1("fa1","sin(x)/x",0,10);
74 
75  // Create an observable
76  RooRealVar x3("x3","x3",0.01,20) ;
77 
78  // Create binding of TF1 object to above observable
79  RooAbsReal* rfa1 = bindFunction(fa1,x3) ;
80 
81  // Print rfa1 definition
82  rfa1->Print() ;
83 
84  // Make plot frame in observable, plot TF1 binding function
85  RooPlot* frame3 = x3.frame(Title("TF1 bound as RooFit function")) ;
86  rfa1->plotOn(frame3) ;
87 
88 
89 
90  TCanvas* c = new TCanvas("rf105_funcbinding","rf105_funcbinding",1200,400) ;
91  c->Divide(3) ;
92  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.6) ; frame1->Draw() ;
93  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
94  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
95 
96 }
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:262
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.
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
RooAbsPdf * bindPdf(const char *name, CFUNCD1D func, RooAbsReal &x)
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:679
double beta_pdf(double x, double a, double b)
Probability density function of the beta distribution.
TArc * a
Definition: textangle.C:12
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 RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Plot (project) PDF on specified frame.
Definition: RooAbsPdf.h:105
double beta(double x, double y)
Calculates the beta function.
RooAbsReal * bindFunction(const char *name, CFUNCD1D func, RooAbsReal &x)
RooCmdArg Title(const char *name)
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition: TMath.cxx:187
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The Canvas class.
Definition: TCanvas.h:31
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1135
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:1702
1-Dim function class
Definition: TF1.h:150
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define gPad
Definition: TVirtualPad.h:284
virtual RooFitResult * fitTo(RooAbsData &data, 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())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1056
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
static const double x3[11]