ROOT logo

From $ROOTSYS/tutorials/roofit/rf101_basics.C

//////////////////////////////////////////////////////////////////////////
//
// 'BASIC FUNCTIONALITY' RooFit tutorial macro #101
// 
// Fitting, plotting, toy data generation on one-dimensional p.d.f
//
// pdf = gauss(x,m,s) 
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "TCanvas.h"
#include "RooPlot.h"
#include "TAxis.h"
using namespace RooFit ;


void rf101_basics()
{
  // S e t u p   m o d e l 
  // ---------------------

  // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
  RooRealVar x("x","x",-10,10) ;
  RooRealVar mean("mean","mean of gaussian",1,-10,10) ;
  RooRealVar sigma("sigma","width of gaussian",1,0.1,10) ;

  // Build gaussian p.d.f in terms of x,mean and sigma
  RooGaussian gauss("gauss","gaussian PDF",x,mean,sigma) ;  

  // Construct plot frame in 'x'
  RooPlot* xframe = x.frame(Title("Gaussian p.d.f.")) ;


  // 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
  // ---------------------------------------------------------------------------

  // Plot gauss in frame (i.e. in x) 
  gauss.plotOn(xframe) ;

  // Change the value of sigma to 3
  sigma.setVal(3) ;

  // Plot gauss in frame (i.e. in x) and draw frame on canvas
  gauss.plotOn(xframe,LineColor(kRed)) ;
  

  // G e n e r a t e   e v e n t s 
  // -----------------------------

  // Generate a dataset of 1000 events in x from gauss
  RooDataSet* data = gauss.generate(x,10000) ;  
  
  // Make a second plot frame in x and draw both the 
  // data and the p.d.f in the frame
  RooPlot* xframe2 = x.frame(Title("Gaussian p.d.f. with data")) ;
  data->plotOn(xframe2) ;
  gauss.plotOn(xframe2) ;
  

  // F i t   m o d e l   t o   d a t a
  // -----------------------------

  // Fit pdf to data
  gauss.fitTo(*data) ;

  // Print values of mean and sigma (that now reflect fitted values and errors)
  mean.Print() ;
  sigma.Print() ;

  // Draw all frames on a canvas
  TCanvas* c = new TCanvas("rf101_basics","rf101_basics",800,400) ;
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.6) ; xframe->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; xframe2->GetYaxis()->SetTitleOffset(1.6) ; xframe2->Draw() ;
  
 
}
 rf101_basics.C:1
 rf101_basics.C:2
 rf101_basics.C:3
 rf101_basics.C:4
 rf101_basics.C:5
 rf101_basics.C:6
 rf101_basics.C:7
 rf101_basics.C:8
 rf101_basics.C:9
 rf101_basics.C:10
 rf101_basics.C:11
 rf101_basics.C:12
 rf101_basics.C:13
 rf101_basics.C:14
 rf101_basics.C:15
 rf101_basics.C:16
 rf101_basics.C:17
 rf101_basics.C:18
 rf101_basics.C:19
 rf101_basics.C:20
 rf101_basics.C:21
 rf101_basics.C:22
 rf101_basics.C:23
 rf101_basics.C:24
 rf101_basics.C:25
 rf101_basics.C:26
 rf101_basics.C:27
 rf101_basics.C:28
 rf101_basics.C:29
 rf101_basics.C:30
 rf101_basics.C:31
 rf101_basics.C:32
 rf101_basics.C:33
 rf101_basics.C:34
 rf101_basics.C:35
 rf101_basics.C:36
 rf101_basics.C:37
 rf101_basics.C:38
 rf101_basics.C:39
 rf101_basics.C:40
 rf101_basics.C:41
 rf101_basics.C:42
 rf101_basics.C:43
 rf101_basics.C:44
 rf101_basics.C:45
 rf101_basics.C:46
 rf101_basics.C:47
 rf101_basics.C:48
 rf101_basics.C:49
 rf101_basics.C:50
 rf101_basics.C:51
 rf101_basics.C:52
 rf101_basics.C:53
 rf101_basics.C:54
 rf101_basics.C:55
 rf101_basics.C:56
 rf101_basics.C:57
 rf101_basics.C:58
 rf101_basics.C:59
 rf101_basics.C:60
 rf101_basics.C:61
 rf101_basics.C:62
 rf101_basics.C:63
 rf101_basics.C:64
 rf101_basics.C:65
 rf101_basics.C:66
 rf101_basics.C:67
 rf101_basics.C:68
 rf101_basics.C:69
 rf101_basics.C:70
 rf101_basics.C:71
 rf101_basics.C:72
 rf101_basics.C:73
 rf101_basics.C:74
 rf101_basics.C:75
 rf101_basics.C:76
 rf101_basics.C:77
 rf101_basics.C:78
 rf101_basics.C:79
 rf101_basics.C:80
 rf101_basics.C:81
 rf101_basics.C:82
 rf101_basics.C:83
 rf101_basics.C:84
 rf101_basics.C:85
 rf101_basics.C:86
 rf101_basics.C:87