Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf407_latextables.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// Data and categories: latex printing of lists and sets of RooArgSets
5 ///
6 /// \macro_output
7 /// \macro_code
8 /// \author 07/2008 - Wouter Verkerke
9 
10 #include "RooRealVar.h"
11 #include "RooDataSet.h"
12 #include "RooGaussian.h"
13 #include "RooConstVar.h"
14 #include "RooChebychev.h"
15 #include "RooAddPdf.h"
16 #include "RooExponential.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "RooPlot.h"
20 using namespace RooFit;
21 
22 void rf407_latextables()
23 {
24  // S e t u p c o m p o s i t e p d f
25  // --------------------------------------
26 
27  // Declare observable x
28  RooRealVar x("x", "x", 0, 10);
29 
30  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
31  RooRealVar mean("mean", "mean of gaussians", 5);
32  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
33  RooRealVar sigma2("sigma2", "width of gaussians", 1);
34  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
35  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
36 
37  // Sum the signal components into a composite signal p.d.f.
38  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
39  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
40 
41  // Build Chebychev polynomial p.d.f.
42  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
43  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
44  RooChebychev bkg1("bkg1", "Background 1", x, RooArgSet(a0, a1));
45 
46  // Build expontential pdf
47  RooRealVar alpha("alpha", "alpha", -1);
48  RooExponential bkg2("bkg2", "Background 2", x, alpha);
49 
50  // Sum the background components into a composite background p.d.f.
51  RooRealVar bkg1frac("sig1frac", "fraction of component 1 in background", 0.2, 0., 1.);
52  RooAddPdf bkg("bkg", "Signal", RooArgList(bkg1, bkg2), sig1frac);
53 
54  // Sum the composite signal and background
55  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
56  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
57 
58  // M a k e l i s t o f p a r a m e t e r s b e f o r e a n d a f t e r f i t
59  // ----------------------------------------------------------------------------------------
60 
61  // Make list of model parameters
62  RooArgSet *params = model.getParameters(x);
63 
64  // Save snapshot of prefit parameters
65  RooArgSet *initParams = (RooArgSet *)params->snapshot();
66 
67  // Do fit to data, to obtain error estimates on parameters
68  RooDataSet *data = model.generate(x, 1000);
69  model.fitTo(*data);
70 
71  // P r i n t l a t ex t a b l e o f p a r a m e t e r s o f p d f
72  // --------------------------------------------------------------------------
73 
74  // Print parameter list in LaTeX for (one column with names, one column with values)
75  params->printLatex();
76 
77  // Print parameter list in LaTeX for (names values|names values)
78  params->printLatex(Columns(2));
79 
80  // Print two parameter lists side by side (name values initvalues)
81  params->printLatex(Sibling(*initParams));
82 
83  // Print two parameter lists side by side (name values initvalues|name values initvalues)
84  params->printLatex(Sibling(*initParams), Columns(2));
85 
86  // Write LaTex table to file
87  params->printLatex(Sibling(*initParams), OutputFile("rf407_latextables.tex"));
88 }