Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf103_interprfuncs.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Basic functionality: interpreted functions and p.d.f.s
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "TCanvas.h"
15 #include "TAxis.h"
16 #include "RooPlot.h"
17 #include "RooFitResult.h"
18 #include "RooGenericPdf.h"
19 #include "RooConstVar.h"
20 
21 using namespace RooFit;
22 
23 void rf103_interprfuncs()
24 {
25  // ----------------------------------------------------
26  // G e n e r i c i n t e r p r e t e d p . d . f .
27  // ====================================================
28 
29  // Declare observable x
30  RooRealVar x("x", "x", -20, 20);
31 
32  // C o n s t r u c t g e n e r i c p d f f r o m i n t e r p r e t e d e x p r e s s i o n
33  // -------------------------------------------------------------------------------------------------
34 
35  // To construct a proper p.d.f, the formula expression is explicitly normalized internally by dividing
36  // it by a numeric integral of the expression over x in the range [-20,20]
37  //
38  RooRealVar alpha("alpha", "alpha", 5, 0.1, 10);
39  RooGenericPdf genpdf("genpdf", "genpdf", "(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))", RooArgSet(x, alpha));
40 
41  // S a m p l e , f i t a n d p l o t g e n e r i c p d f
42  // ---------------------------------------------------------------
43 
44  // Generate a toy dataset from the interpreted p.d.f
45  RooDataSet *data = genpdf.generate(x, 10000);
46 
47  // Fit the interpreted p.d.f to the generated data
48  genpdf.fitTo(*data);
49 
50  // Make a plot of the data and the p.d.f overlaid
51  RooPlot *xframe = x.frame(Title("Interpreted expression pdf"));
52  data->plotOn(xframe);
53  genpdf.plotOn(xframe);
54 
55  // -----------------------------------------------------------------------------------------------------------
56  // S t a n d a r d p . d . f a d j u s t w i t h i n t e r p r e t e d h e l p e r f u n c t i o n
57  // ==========================================================================================================
58  // Make a gauss(x,sqrt(mean2),sigma) from a standard RooGaussian
59 
60  // C o n s t r u c t s t a n d a r d p d f w i t h f o r m u l a r e p l a c i n g p a r a m e t e r
61  // ------------------------------------------------------------------------------------------------------------
62 
63  // Construct parameter mean2 and sigma
64  RooRealVar mean2("mean2", "mean^2", 10, 0, 200);
65  RooRealVar sigma("sigma", "sigma", 3, 0.1, 10);
66 
67  // Construct interpreted function mean = sqrt(mean^2)
68  RooFormulaVar mean("mean", "mean", "sqrt(mean2)", mean2);
69 
70  // Construct a gaussian g2(x,sqrt(mean2),sigma) ;
71  RooGaussian g2("g2", "h2", x, mean, sigma);
72 
73  // G e n e r a t e t o y d a t a
74  // ---------------------------------
75 
76  // Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian dataset with mean 10 and width 3
77  RooGaussian g1("g1", "g1", x, RooConst(10), RooConst(3));
78  RooDataSet *data2 = g1.generate(x, 1000);
79 
80  // F i t a n d p l o t t a i l o r e d s t a n d a r d p d f
81  // -------------------------------------------------------------------
82 
83  // Fit g2 to data from g1
84  RooFitResult *r = g2.fitTo(*data2, Save());
85  r->Print();
86 
87  // Plot data on frame and overlay projection of g2
88  RooPlot *xframe2 = x.frame(Title("Tailored Gaussian pdf"));
89  data2->plotOn(xframe2);
90  g2.plotOn(xframe2);
91 
92  // Draw all frames on a canvas
93  TCanvas *c = new TCanvas("rf103_interprfuncs", "rf103_interprfuncs", 800, 400);
94  c->Divide(2);
95  c->cd(1);
96  gPad->SetLeftMargin(0.15);
97  xframe->GetYaxis()->SetTitleOffset(1.4);
98  xframe->Draw();
99  c->cd(2);
100  gPad->SetLeftMargin(0.15);
101  xframe2->GetYaxis()->SetTitleOffset(1.4);
102  xframe2->Draw();
103 }