Data and categories: latex printing of lists and sets of RooArgSets
import ROOT
x = ROOT.RooRealVar("x", "x", 0, 10)
mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
bkg1 = ROOT.RooChebychev("bkg1", "Background 1", x, [a0, a1])
alpha = ROOT.RooRealVar("alpha", "alpha", -1)
bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
bkg1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [sig1frac])
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
params = model.getParameters({x})
initParams = params.snapshot()
data = model.generate({x}, 1000)
model.fitTo(data, PrintLevel=-1)
params.printLatex()
params.printLatex(Columns=2)
params.printLatex(Sibling=initParams)
params.printLatex(Sibling=initParams, Columns=2)
params.printLatex(Sibling=initParams, OutputFile="rf407_latextables.tex")
[#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.
[#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
\begin{tabular}{lc}
$\verb+a0+ $ & $ 0.5\pm 0.2$\\
$\verb+a1+ $ & $ 0.3\pm 0.1$\\
$\verb+alpha+ $ & $ -1.00$\\
$\verb+bkgfrac+ $ & $ 0.46\pm 0.03$\\
$\verb+mean+ $ & $ 5$\\
$\verb+sig1frac+ $ & $ 0.79\pm 0.05$\\
$\verb+sigma1+ $ & $ 0.5$\\
$\verb+sigma2+ $ & $ 1$\\
\end{tabular}
\begin{tabular}{lc|lc}
$\verb+a0+ $ & $ 0.5\pm 0.2$ & $\verb+mean+ $ & $ 5$\\
$\verb+a1+ $ & $ 0.3\pm 0.1$ & $\verb+sig1frac+ $ & $ 0.79\pm 0.05$\\
$\verb+alpha+ $ & $ -1.00$ & $\verb+sigma1+ $ & $ 0.5$\\
$\verb+bkgfrac+ $ & $ 0.46\pm 0.03$ & $\verb+sigma2+ $ & $ 1$\\
\end{tabular}
\begin{tabular}{lcc}
$\verb+a0+ $ & $ 0.5\pm 0.2$ & $ 0.5$\\
$\verb+a1+ $ & $ 0.3\pm 0.1$ & $ 0$\\
$\verb+alpha+ $ & $ -1.00$ & $-1.00$\\
$\verb+bkgfrac+ $ & $ 0.46\pm 0.03$ & $ 0.5$\\
$\verb+mean+ $ & $ 5$ & $ 5$\\
$\verb+sig1frac+ $ & $ 0.79\pm 0.05$ & $ 0.8$\\
$\verb+sigma1+ $ & $ 0.5$ & $ 0.5$\\
$\verb+sigma2+ $ & $ 1$ & $ 1$\\
\end{tabular}
\begin{tabular}{lcc|lcc}
$\verb+a0+ $ & $ 0.5\pm 0.2$ & $ 0.5$ & $\verb+mean+ $ & $ 5$ & $ 5$\\
$\verb+a1+ $ & $ 0.3\pm 0.1$ & $ 0$ & $\verb+sig1frac+ $ & $ 0.79\pm 0.05$ & $ 0.8$\\
$\verb+alpha+ $ & $ -1.00$ & $-1.00$ & $\verb+sigma1+ $ & $ 0.5$ & $ 0.5$\\
$\verb+bkgfrac+ $ & $ 0.46\pm 0.03$ & $ 0.5$ & $\verb+sigma2+ $ & $ 1$ & $ 1$\\
\end{tabular}
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf407_latextables.py.