Data and categories: latex printing of lists and sets of RooArgSets
Author: Wouter Verkerke
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, April 17, 2024 at 11:18 AM.
%%cpp -d
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooChebychev.h"
#include "RooAddPdf.h"
#include "RooExponential.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;
Declare observable x
RooRealVar x("x", "x", 0, 10);
Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
RooRealVar mean("mean", "mean of gaussians", 5);
RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
RooRealVar sigma2("sigma2", "width of gaussians", 1);
RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-inf, inf] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range. [#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-inf, inf] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
Sum the signal components into a composite signal pdf
RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
Build Chebychev polynomial pdf
RooRealVar a0("a0", "a0", 0.5, 0., 1.);
RooRealVar a1("a1", "a1", 0.2, 0., 1.);
RooChebychev bkg1("bkg1", "Background 1", x, RooArgSet(a0, a1));
Build expontential pdf
RooRealVar alpha("alpha", "alpha", -1);
RooExponential bkg2("bkg2", "Background 2", x, alpha);
Sum the background components into a composite background pdf
RooRealVar bkg1frac("sig1frac", "fraction of component 1 in background", 0.2, 0., 1.);
RooAddPdf bkg("bkg", "Signal", RooArgList(bkg1, bkg2), sig1frac);
Sum the composite signal and background
RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
Make list of model parameters
std::unique_ptr<RooArgSet> params{model.getParameters(x)};
Save snapshot of prefit parameters
std::unique_ptr<RooArgSet> initParams{static_cast<RooArgSet *>(params->snapshot())};
Do fit to data, to obtain error estimates on parameters
std::unique_ptr<RooDataSet> data{model.generate(x, 1000)};
model.fitTo(*data, PrintLevel(-1));
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data [#1] INFO:Fitting -- using CPU computation library compiled with -mavx2 [#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
input_line_57:2:2: warning: 'data' shadows a declaration with the same name in the 'std' namespace; use '::data' to reference this declaration std::unique_ptr<RooDataSet> data{model.generate(x, 1000)}; ^
Print parameter list in LaTeX for (one column with names, one column with values)
params->printLatex();
\begin{tabular}{lc} $\verb+a0+ $ & $ 0.6\pm 0.2$\\ $\verb+a1+ $ & $ 0.2\pm 0.2$\\ $\verb+alpha+ $ & $ -1.00$\\ $\verb+bkgfrac+ $ & $ 0.45\pm 0.03$\\ $\verb+mean+ $ & $ 5$\\ $\verb+sig1frac+ $ & $ 0.71\pm 0.06$\\ $\verb+sigma1+ $ & $ 0.5$\\ $\verb+sigma2+ $ & $ 1$\\ \end{tabular}
Print parameter list in LaTeX for (names values|names values)
params->printLatex(Columns(2));
\begin{tabular}{lc|lc} $\verb+a0+ $ & $ 0.6\pm 0.2$ & $\verb+mean+ $ & $ 5$\\ $\verb+a1+ $ & $ 0.2\pm 0.2$ & $\verb+sig1frac+ $ & $ 0.71\pm 0.06$\\ $\verb+alpha+ $ & $ -1.00$ & $\verb+sigma1+ $ & $ 0.5$\\ $\verb+bkgfrac+ $ & $ 0.45\pm 0.03$ & $\verb+sigma2+ $ & $ 1$\\ \end{tabular}
Print two parameter lists side by side (name values initvalues)
params->printLatex(Sibling(*initParams));
\begin{tabular}{lcc} $\verb+a0+ $ & $ 0.6\pm 0.2$ & $ 0.5$\\ $\verb+a1+ $ & $ 0.2\pm 0.2$ & $ 0.2$\\ $\verb+alpha+ $ & $ -1.00$ & $-1.00$\\ $\verb+bkgfrac+ $ & $ 0.45\pm 0.03$ & $ 0.5$\\ $\verb+mean+ $ & $ 5$ & $ 5$\\ $\verb+sig1frac+ $ & $ 0.71\pm 0.06$ & $ 0.8$\\ $\verb+sigma1+ $ & $ 0.5$ & $ 0.5$\\ $\verb+sigma2+ $ & $ 1$ & $ 1$\\ \end{tabular}
Print two parameter lists side by side (name values initvalues|name values initvalues)
params->printLatex(Sibling(*initParams), Columns(2));
\begin{tabular}{lcc|lcc} $\verb+a0+ $ & $ 0.6\pm 0.2$ & $ 0.5$ & $\verb+mean+ $ & $ 5$ & $ 5$\\ $\verb+a1+ $ & $ 0.2\pm 0.2$ & $ 0.2$ & $\verb+sig1frac+ $ & $ 0.71\pm 0.06$ & $ 0.8$\\ $\verb+alpha+ $ & $ -1.00$ & $-1.00$ & $\verb+sigma1+ $ & $ 0.5$ & $ 0.5$\\ $\verb+bkgfrac+ $ & $ 0.45\pm 0.03$ & $ 0.5$ & $\verb+sigma2+ $ & $ 1$ & $ 1$\\ \end{tabular}
Write LaTex table to file
params->printLatex(Sibling(*initParams), OutputFile("rf407_latextables.tex"));