Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf101_basics.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Basic functionality: fitting, plotting, toy data generation on one-dimensional p.d.f
5 ///
6 /// pdf = gauss(x,m,s)
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "TCanvas.h"
17 #include "RooPlot.h"
18 #include "TAxis.h"
19 using namespace RooFit;
20 
21 void rf101_basics()
22 {
23  // S e t u p m o d e l
24  // ---------------------
25 
26  // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
27  RooRealVar x("x", "x", -10, 10);
28  RooRealVar mean("mean", "mean of gaussian", 1, -10, 10);
29  RooRealVar sigma("sigma", "width of gaussian", 1, 0.1, 10);
30 
31  // Build gaussian p.d.f in terms of x,mean and sigma
32  RooGaussian gauss("gauss", "gaussian PDF", x, mean, sigma);
33 
34  // Construct plot frame in 'x'
35  RooPlot *xframe = x.frame(Title("Gaussian p.d.f."));
36 
37  // 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
38  // ---------------------------------------------------------------------------
39 
40  // Plot gauss in frame (i.e. in x)
41  gauss.plotOn(xframe);
42 
43  // Change the value of sigma to 3
44  sigma.setVal(3);
45 
46  // Plot gauss in frame (i.e. in x) and draw frame on canvas
47  gauss.plotOn(xframe, LineColor(kRed));
48 
49  // G e n e r a t e e v e n t s
50  // -----------------------------
51 
52  // Generate a dataset of 1000 events in x from gauss
53  RooDataSet *data = gauss.generate(x, 10000);
54 
55  // Make a second plot frame in x and draw both the
56  // data and the p.d.f in the frame
57  RooPlot *xframe2 = x.frame(Title("Gaussian p.d.f. with data"));
58  data->plotOn(xframe2);
59  gauss.plotOn(xframe2);
60 
61  // F i t m o d e l t o d a t a
62  // -----------------------------
63 
64  // Fit pdf to data
65  gauss.fitTo(*data);
66 
67  // Print values of mean and sigma (that now reflect fitted values and errors)
68  mean.Print();
69  sigma.Print();
70 
71  // Draw all frames on a canvas
72  TCanvas *c = new TCanvas("rf101_basics", "rf101_basics", 800, 400);
73  c->Divide(2);
74  c->cd(1);
75  gPad->SetLeftMargin(0.15);
76  xframe->GetYaxis()->SetTitleOffset(1.6);
77  xframe->Draw();
78  c->cd(2);
79  gPad->SetLeftMargin(0.15);
80  xframe2->GetYaxis()->SetTitleOffset(1.6);
81  xframe2->Draw();
82 }