Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
rf407_ComputationalGraphVisualization.py File Reference

Detailed Description

View in nbviewer Open in SWAN
Data and categories: Visualing computational graph model before fitting, and latex printing of lists and sets of RooArgSets after fitting

import ROOT
# Setup composite pdf
# --------------------------------------
# Declare observable x
x = ROOT.RooRealVar("x", "x", 0, 10)
# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
# their parameters
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)
# Sum the signal components into a composite signal pdf
sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
# Build Chebychev polynomial pdf
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])
# Build expontential pdf
alpha = ROOT.RooRealVar("alpha", "alpha", -1)
bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
# Sum the background components into a composite background pdf
bkg1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [sig1frac])
# Sum the composite signal and background
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
# Print composite tree in ASCII
# -----------------------------------------------------------
# Print tree to stdout
# Print tree to file
model.printCompactTree("", "rf206_asciitree.txt")
# Draw composite tree graphically
# -------------------------------------------------------------
# Print GraphViz DOT file with representation of tree
model.graphVizTree("rf206_model.dot")
# Make list of parameters before and after fit
# ----------------------------------------------------------------------------------------
# Make list of model parameters
params = model.getParameters({x})
# Save snapshot of prefit parameters
initParams = params.snapshot()
# Do fit to data, obtain error estimates on parameters
data = model.generate({x}, 1000)
model.fitTo(data, PrintLevel=-1)
# Print LateX table of parameters of pdf
# --------------------------------------------------------------------------
# Print parameter list in LaTeX for (one column with names, column with
# values)
# Print parameter list in LaTeX for (names values|names values)
# Print two parameter lists side by side (name values initvalues)
params.printLatex(Sibling=initParams)
# Print two parameter lists side by side (name values initvalues|name
# values initvalues)
params.printLatex(Sibling=initParams, Columns=2)
# Write LaTex table to file
params.printLatex(Sibling=initParams, OutputFile="rf407_latextables.tex")
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
[#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.
0x55ac3cf5fe50 RooAddPdf::model = 0.900674/1 [Auto,Clean]
0x55ac3ce1fd90/V- RooAddPdf::bkg = 0.801348/1 [Auto,Clean]
0x55ac3ce19070/V- RooChebychev::bkg1 = 1 [Auto,Dirty]
0x55ac3c5bc330/V- RooRealVar::x = 5
0x55ac3ca88ef0/V- RooRealVar::a0 = 0.5
0x55ac3caa82d0/V- RooRealVar::a1 = 0
0x55ac3ca29770/V- RooRealVar::sig1frac = 0.8
0x55ac3cdef020/V- RooExponential::bkg2 = 0.00673795 [Auto,Dirty]
0x55ac3c5bc330/V- RooRealVar::x = 5
0x55ac3c9db080/V- RooRealVar::alpha = -1
0x55ac38585e20/V- RooRealVar::bkgfrac = 0.5
0x55ac3a10ac00/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x55ac3ca51a60/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x55ac3c5bc330/V- RooRealVar::x = 5
0x55ac3c70bb20/V- RooRealVar::mean = 5
0x55ac3c50c380/V- RooRealVar::sigma1 = 0.5
0x55ac3ca29770/V- RooRealVar::sig1frac = 0.8
0x55ac3c936eb0/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x55ac3c5bc330/V- RooRealVar::x = 5
0x55ac3c70bb20/V- RooRealVar::mean = 5
0x55ac387ce8b0/V- RooRealVar::sigma2 = 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 -mavx512
[#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_ComputationalGraphVisualization.py.