Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf206_treevistools.py File Reference

Detailed Description

View in nbviewer Open in SWAN
Addition and convolution: tools for visualization of ROOT.RooAbsArg expression trees

import ROOT
# Set up 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("bkg1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [bkg1frac])
# 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
model.Print("t")
# 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 graphic output file with one of the GraphViz tools
# (freely available from www.graphviz.org)
#
# 'Top-to-bottom graph'
# unix> dot -Tgif -o rf207_model_dot.gif rf207_model.dot
#
# 'Spring-model graph'
# unix> fdp -Tgif -o rf207_model_fdp.gif rf207_model.dot
[#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.
0x89f2d90 RooAddPdf::model = 0.602695/1 [Auto,Clean]
0x89f1c10/V- RooAddPdf::bkg = 0.20539/1 [Auto,Clean]
0x8487580/V- RooChebychev::bkg1 = 1 [Auto,Dirty]
0x8253930/V- RooRealVar::x = 5
0x88e0ac0/V- RooRealVar::a0 = 0.5
0x5cfc070/V- RooRealVar::a1 = 0
0x89ec6d0/V- RooRealVar::bkg1frac = 0.2
0x89f9750/V- RooExponential::bkg2 = 0.00673795 [Auto,Dirty]
0x8253930/V- RooRealVar::x = 5
0x5cfa940/V- RooRealVar::alpha = -1
0x88ecb80/V- RooRealVar::bkgfrac = 0.5
0x8996f60/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x84999d0/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x8253930/V- RooRealVar::x = 5
0x80f9110/V- RooRealVar::mean = 5
0x7ddb290/V- RooRealVar::sigma1 = 0.5
0x68fc3c0/V- RooRealVar::sig1frac = 0.8
0x84ab360/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x8253930/V- RooRealVar::x = 5
0x80f9110/V- RooRealVar::mean = 5
0x7d23ae0/V- RooRealVar::sigma2 = 1
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf206_treevistools.py.