Logo ROOT   6.10/09
Reference Guide
rf105_funcbinding.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #105
4 //
5 // Demonstration of binding ROOT Math functions as RooFit functions
6 // and pdfs
7 //
8 // 07/2008 - Wouter Verkerke
9 //
10 /////////////////////////////////////////////////////////////////////////
11 
12 #ifndef __CINT__
13 #include "RooGlobalFunc.h"
14 #endif
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "TCanvas.h"
19 #include "RooPlot.h"
20 #include "TMath.h"
21 #include "TF1.h"
22 #include "Math/DistFunc.h"
23 #include "RooCFunction1Binding.h"
24 #include "RooCFunction3Binding.h"
25 #include "RooTFnBinding.h"
26 
27 using namespace RooFit ;
28 
29 class TestBasic105 : public RooFitTestUnit
30 {
31 public:
32  TestBasic105(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("C++ function binding operator p.d.f",refFile,writeRef,verbose) {} ;
33  Bool_t testCode() {
34 
35  // B i n d T M a t h : : E r f C f u n c t i o n
36  // ---------------------------------------------------
37 
38  // Bind one-dimensional TMath::Erf function as RooAbsReal function
39  RooRealVar x("x","x",-3,3) ;
41 
42  // Plot erf on frame
43  RooPlot* frame1 = x.frame(Title("TMath::Erf bound as RooFit function")) ;
44  erf->plotOn(frame1) ;
45 
46 
47  // 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
48  // -----------------------------------------------------------------------
49 
50  // Bind pdf ROOT::Math::Beta with three variables as RooAbsPdf function
51  RooRealVar x2("x2","x2",0,1) ;
52  RooRealVar a("a","a",5,0,10) ;
53  RooRealVar b("b","b",2,0,10) ;
55 
56  // Generate some events and fit
57  RooDataSet* data = beta->generate(x2,10000) ;
58  beta->fitTo(*data) ;
59 
60  // Plot data and pdf on frame
61  RooPlot* frame2 = x2.frame(Title("ROOT::Math::Beta bound as RooFit pdf")) ;
62  data->plotOn(frame2) ;
63  beta->plotOn(frame2) ;
64 
65 
66 
67  // 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
68  // ---------------------------------------------------------------
69 
70  // Create a ROOT TF1 function
71  TF1 *fa1 = new TF1("fa1","sin(x)/x",0,10);
72 
73  // Create an observable
74  RooRealVar x3("x3","x3",0.01,20) ;
75 
76  // Create binding of TF1 object to above observable
77  RooAbsReal* rfa1 = bindFunction(fa1,x3) ;
78 
79  // Make plot frame in observable, plot TF1 binding function
80  RooPlot* frame3 = x3.frame(Title("TF1 bound as RooFit function")) ;
81  rfa1->plotOn(frame3) ;
82 
83 
84  regPlot(frame1,"rf105_plot1") ;
85  regPlot(frame2,"rf105_plot2") ;
86  regPlot(frame3,"rf105_plot3") ;
87 
88  delete erf ;
89  delete beta ;
90  delete fa1 ;
91  delete rfa1 ;
92  delete data ;
93 
94  return kTRUE ;
95  }
96 } ;
double erf(double x)
Error function encountered in integrating the normal distribution.
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.
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)
double beta_pdf(double x, double a, double b)
Probability density function of the beta distribution.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
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
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
bool verbose
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
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
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
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
const Bool_t kTRUE
Definition: RtypesCore.h:91
static const double x3[11]