Logo ROOT   6.10/09
Reference Guide
rf101_basics.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #101
4 //
5 // Fitting, plotting, toy data generation on one-dimensional p.d.f
6 //
7 // pdf = gauss(x,m,s)
8 //
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooDataSet.h"
19 #include "RooGaussian.h"
20 #include "TCanvas.h"
21 #include "RooPlot.h"
22 using namespace RooFit ;
23 
24 
25 // Elementary operations on a gaussian PDF
26 class TestBasic101 : public RooFitTestUnit
27 {
28 public:
29  TestBasic101(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fitting,plotting & event generation of basic p.d.f",refFile,writeRef,verbose) {} ;
30  Bool_t testCode() {
31 
32  // S e t u p m o d e l
33  // ---------------------
34 
35  // Declare variables x,mean,sigma with associated name, title, initial value and allowed range
36  RooRealVar x("x","x",-10,10) ;
37  RooRealVar mean("mean","mean of gaussian",1,-10,10) ;
38  RooRealVar sigma("sigma","width of gaussian",1,0.1,10) ;
39 
40  // Build gaussian p.d.f in terms of x,mean and sigma
41  RooGaussian gauss("gauss","gaussian PDF",x,mean,sigma) ;
42 
43  // Construct plot frame in 'x'
44  RooPlot* xframe = x.frame(Title("Gaussian p.d.f.")) ;
45 
46 
47  // 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
48  // ---------------------------------------------------------------------------
49 
50  // Plot gauss in frame (i.e. in x)
51  gauss.plotOn(xframe) ;
52 
53  // Change the value of sigma to 3
54  sigma.setVal(3) ;
55 
56  // Plot gauss in frame (i.e. in x) and draw frame on canvas
57  gauss.plotOn(xframe,LineColor(kRed),Name("another")) ;
58 
59 
60  // G e n e r a t e e v e n t s
61  // -----------------------------
62 
63  // Generate a dataset of 1000 events in x from gauss
64  RooDataSet* data = gauss.generate(x,10000) ;
65 
66  // Make a second plot frame in x and draw both the
67  // data and the p.d.f in the frame
68  RooPlot* xframe2 = x.frame(Title("Gaussian p.d.f. with data")) ;
69  data->plotOn(xframe2) ;
70  gauss.plotOn(xframe2) ;
71 
72 
73  // F i t m o d e l t o d a t a
74  // -----------------------------
75 
76  // Fit pdf to data
77  gauss.fitTo(*data) ;
78 
79 
80  // --- Post processing for stressRooFit ---
81  regPlot(xframe ,"rf101_plot1") ;
82  regPlot(xframe2,"rf101_plot2") ;
83 
84  delete data ;
85 
86  return kTRUE ;
87  }
88 } ;
89 
90 
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
const Double_t sigma
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooCmdArg Name(const char *name)
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
const Bool_t kTRUE
Definition: RtypesCore.h:91