rf101_basics.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #101
00004 // 
00005 // Fitting, plotting, toy data generation on one-dimensional p.d.f
00006 //
00007 // pdf = gauss(x,m,s) 
00008 //
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 // 
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "TCanvas.h"
00021 #include "RooPlot.h"
00022 using namespace RooFit ;
00023 
00024 
00025 // Elementary operations on a gaussian PDF
00026 class TestBasic101 : public RooFitTestUnit
00027 {
00028 public: 
00029   TestBasic101(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fitting,plotting & event generation of basic p.d.f",refFile,writeRef,verbose) {} ;
00030   Bool_t testCode() {
00031 
00032     // S e t u p   m o d e l 
00033     // ---------------------
00034     
00035     // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
00036     RooRealVar x("x","x",-10,10) ;
00037     RooRealVar mean("mean","mean of gaussian",1,-10,10) ;
00038     RooRealVar sigma("sigma","width of gaussian",1,0.1,10) ;
00039     
00040     // Build gaussian p.d.f in terms of x,mean and sigma
00041     RooGaussian gauss("gauss","gaussian PDF",x,mean,sigma) ;  
00042     
00043     // Construct plot frame in 'x'
00044     RooPlot* xframe = x.frame(Title("Gaussian p.d.f.")) ;
00045     
00046     
00047     // P l o t   m o d e l   a n d   c h a n g e   p a r a m e t e r   v a l u e s
00048     // ---------------------------------------------------------------------------
00049     
00050     // Plot gauss in frame (i.e. in x) 
00051     gauss.plotOn(xframe) ;
00052     
00053     // Change the value of sigma to 3
00054     sigma.setVal(3) ;
00055     
00056     // Plot gauss in frame (i.e. in x) and draw frame on canvas
00057     gauss.plotOn(xframe,LineColor(kRed),Name("another")) ;
00058     
00059     
00060     // G e n e r a t e   e v e n t s 
00061     // -----------------------------
00062     
00063     // Generate a dataset of 1000 events in x from gauss
00064     RooDataSet* data = gauss.generate(x,10000) ;  
00065     
00066     // Make a second plot frame in x and draw both the 
00067     // data and the p.d.f in the frame
00068     RooPlot* xframe2 = x.frame(Title("Gaussian p.d.f. with data")) ;
00069     data->plotOn(xframe2) ;
00070     gauss.plotOn(xframe2) ;
00071     
00072     
00073     // F i t   m o d e l   t o   d a t a
00074     // -----------------------------
00075     
00076     // Fit pdf to data
00077     gauss.fitTo(*data) ;
00078         
00079     
00080     // --- Post processing for stressRooFit ---
00081     regPlot(xframe ,"rf101_plot1") ;
00082     regPlot(xframe2,"rf101_plot2") ;
00083     
00084     delete data ;
00085     
00086     return kTRUE ;
00087   }
00088 } ;
00089 
00090 

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