rf105_funcbinding.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #105
00004 // 
00005 //  Demonstration of binding ROOT Math functions as RooFit functions
00006 //  and pdfs
00007 //
00008 // 07/2008 - Wouter Verkerke 
00009 // 
00010 /////////////////////////////////////////////////////////////////////////
00011 
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooGaussian.h"
00018 #include "TCanvas.h"
00019 #include "RooPlot.h"
00020 #include "TMath.h"
00021 #include "TF1.h"
00022 #include "Math/DistFunc.h"
00023 #include "RooCFunction1Binding.h" 
00024 #include "RooCFunction3Binding.h" 
00025 #include "RooTFnBinding.h" 
00026 
00027 using namespace RooFit ;
00028 
00029 class TestBasic105 : public RooFitTestUnit
00030 {
00031 public: 
00032   TestBasic105(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("C++ function binding operator p.d.f",refFile,writeRef,verbose) {} ;
00033   Bool_t testCode() {
00034 
00035     // B i n d   T M a t h : : E r f   C   f u n c t i o n
00036     // ---------------------------------------------------
00037     
00038     // Bind one-dimensional TMath::Erf function as RooAbsReal function
00039     RooRealVar x("x","x",-3,3) ;
00040     RooAbsReal* erf = bindFunction("erf",TMath::Erf,x) ;
00041     
00042     // Plot erf on frame 
00043     RooPlot* frame1 = x.frame(Title("TMath::Erf bound as RooFit function")) ;
00044     erf->plotOn(frame1) ;
00045     
00046     
00047     // 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
00048     // -----------------------------------------------------------------------
00049     
00050     // Bind pdf ROOT::Math::Beta with three variables as RooAbsPdf function
00051     RooRealVar x2("x2","x2",0,1) ;
00052     RooRealVar a("a","a",5,0,10) ;
00053     RooRealVar b("b","b",2,0,10) ;
00054     RooAbsPdf* beta = bindPdf("beta",ROOT::Math::beta_pdf,x2,a,b) ;
00055     
00056     // Generate some events and fit
00057     RooDataSet* data = beta->generate(x2,10000) ;
00058     beta->fitTo(*data) ;
00059     
00060     // Plot data and pdf on frame
00061     RooPlot* frame2 = x2.frame(Title("ROOT::Math::Beta bound as RooFit pdf")) ;
00062     data->plotOn(frame2) ;
00063     beta->plotOn(frame2) ;
00064     
00065     
00066     
00067     // 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
00068     // ---------------------------------------------------------------
00069     
00070     // Create a ROOT TF1 function
00071     TF1 *fa1 = new TF1("fa1","sin(x)/x",0,10);
00072     
00073     // Create an observable 
00074     RooRealVar x3("x3","x3",0.01,20) ;
00075     
00076     // Create binding of TF1 object to above observable
00077     RooAbsReal* rfa1 = bindFunction(fa1,x3) ;
00078     
00079     // Make plot frame in observable, plot TF1 binding function
00080     RooPlot* frame3 = x3.frame(Title("TF1 bound as RooFit function")) ;
00081     rfa1->plotOn(frame3) ;
00082       
00083       
00084     regPlot(frame1,"rf105_plot1") ;
00085     regPlot(frame2,"rf105_plot2") ;
00086     regPlot(frame3,"rf105_plot3") ;
00087     
00088     delete erf ;
00089     delete beta ;
00090     delete fa1 ;
00091     delete rfa1 ;
00092     delete data ;
00093     
00094     return kTRUE ;
00095   }  
00096 } ;

Generated on Tue Jul 5 15:07:45 2011 for ROOT_528-00b_version by  doxygen 1.5.1