rf101_basics.C

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 #include "TAxis.h"
00023 using namespace RooFit ;
00024 
00025 
00026 void rf101_basics()
00027 {
00028   // S e t u p   m o d e l 
00029   // ---------------------
00030 
00031   // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
00032   RooRealVar x("x","x",-10,10) ;
00033   RooRealVar mean("mean","mean of gaussian",1,-10,10) ;
00034   RooRealVar sigma("sigma","width of gaussian",1,0.1,10) ;
00035 
00036   // Build gaussian p.d.f in terms of x,mean and sigma
00037   RooGaussian gauss("gauss","gaussian PDF",x,mean,sigma) ;  
00038 
00039   // Construct plot frame in 'x'
00040   RooPlot* xframe = x.frame(Title("Gaussian p.d.f.")) ;
00041 
00042 
00043   // 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
00044   // ---------------------------------------------------------------------------
00045 
00046   // Plot gauss in frame (i.e. in x) 
00047   gauss.plotOn(xframe) ;
00048 
00049   // Change the value of sigma to 3
00050   sigma.setVal(3) ;
00051 
00052   // Plot gauss in frame (i.e. in x) and draw frame on canvas
00053   gauss.plotOn(xframe,LineColor(kRed)) ;
00054   
00055 
00056   // G e n e r a t e   e v e n t s 
00057   // -----------------------------
00058 
00059   // Generate a dataset of 1000 events in x from gauss
00060   RooDataSet* data = gauss.generate(x,10000) ;  
00061   
00062   // Make a second plot frame in x and draw both the 
00063   // data and the p.d.f in the frame
00064   RooPlot* xframe2 = x.frame(Title("Gaussian p.d.f. with data")) ;
00065   data->plotOn(xframe2) ;
00066   gauss.plotOn(xframe2) ;
00067   
00068 
00069   // F i t   m o d e l   t o   d a t a
00070   // -----------------------------
00071 
00072   // Fit pdf to data
00073   gauss.fitTo(*data) ;
00074 
00075   // Print values of mean and sigma (that now reflect fitted values and errors)
00076   mean.Print() ;
00077   sigma.Print() ;
00078 
00079   // Draw all frames on a canvas
00080   TCanvas* c = new TCanvas("rf101_basics","rf101_basics",800,400) ;
00081   c->Divide(2) ;
00082   c->cd(1) ; gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.6) ; xframe->Draw() ;
00083   c->cd(2) ; gPad->SetLeftMargin(0.15) ; xframe2->GetYaxis()->SetTitleOffset(1.6) ; xframe2->Draw() ;
00084   
00085  
00086 }

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