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