Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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///
9/// \date July 2008
10/// \author Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooDataSet.h"
14#include "RooGaussian.h"
15#include "RooConstVar.h"
16#include "RooChebychev.h"
17#include "RooAddPdf.h"
18#include "RooExponential.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
25{
26 // S e t u p c o m p o s i t e p d f
27 // --------------------------------------
28
29 // Declare observable x
30 RooRealVar x("x", "x", 0, 10);
31
32 // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
33 RooRealVar mean("mean", "mean of gaussians", 5);
34 RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
35 RooRealVar sigma2("sigma2", "width of gaussians", 1);
36 RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
37 RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
38
39 // Sum the signal components into a composite signal pdf
40 RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
41 RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
42
43 // Build Chebychev polynomial pdf
44 RooRealVar a0("a0", "a0", 0.5, 0., 1.);
45 RooRealVar a1("a1", "a1", 0.2, 0., 1.);
46 RooChebychev bkg1("bkg1", "Background 1", x, RooArgSet(a0, a1));
47
48 // Build expontential pdf
49 RooRealVar alpha("alpha", "alpha", -1);
50 RooExponential bkg2("bkg2", "Background 2", x, alpha);
51
52 // Sum the background components into a composite background pdf
53 RooRealVar bkg1frac("sig1frac", "fraction of component 1 in background", 0.2, 0., 1.);
54 RooAddPdf bkg("bkg", "Signal", RooArgList(bkg1, bkg2), sig1frac);
55
56 // Sum the composite signal and background
57 RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
58 RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
59
60 // 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
61 // ----------------------------------------------------------------------------------------
62
63 // Make list of model parameters
64 RooArgSet *params = model.getParameters(x);
65
66 // Save snapshot of prefit parameters
67 RooArgSet *initParams = (RooArgSet *)params->snapshot();
68
69 // Do fit to data, to obtain error estimates on parameters
70 RooDataSet *data = model.generate(x, 1000);
71 model.fitTo(*data);
72
73 // 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
74 // --------------------------------------------------------------------------
75
76 // Print parameter list in LaTeX for (one column with names, one column with values)
77 params->printLatex();
78
79 // Print parameter list in LaTeX for (names values|names values)
80 params->printLatex(Columns(2));
81
82 // Print two parameter lists side by side (name values initvalues)
83 params->printLatex(Sibling(*initParams));
84
85 // Print two parameter lists side by side (name values initvalues|name values initvalues)
86 params->printLatex(Sibling(*initParams), Columns(2));
87
88 // Write LaTex table to file
89 params->printLatex(Sibling(*initParams), OutputFile("rf407_latextables.tex"));
90}
void printLatex(const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg()) const
Output content of collection as LaTex table.
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:118
Chebychev polynomial p.d.f.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
Exponential PDF.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
RooCmdArg Columns(Int_t ncol)
RooCmdArg Sibling(const RooAbsCollection &sibling)
RooCmdArg OutputFile(const char *fileName)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...