Addition and convolution: tools for visualization of ROOT.RooAbsArg expression trees
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("bkg1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [bkg1frac])
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
model.Print("t")
model.printCompactTree("", "rf206_asciitree.txt")
model.graphVizTree("rf206_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.
0x7d68ce0 RooAddPdf::model = 0.602695/1 [Auto,Clean]
0x51b2600/V- RooAddPdf::bkg = 0.20539/1 [Auto,Clean]
0x51980c0/V- RooChebychev::bkg1 = 1 [Auto,Dirty]
0x73834a0/V- RooRealVar::x = 5
0x78fc490/V- RooRealVar::a0 = 0.5
0x790bad0/V- RooRealVar::a1 = 0
0x78049b0/V- RooRealVar::bkg1frac = 0.2
0x51b2db0/V- RooExponential::bkg2 = 0.00673795 [Auto,Dirty]
0x73834a0/V- RooRealVar::x = 5
0x7942380/V- RooRealVar::alpha = -1
0x32f07d0/V- RooRealVar::bkgfrac = 0.5
0x7c77c70/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x791d640/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x73834a0/V- RooRealVar::x = 5
0x72df100/V- RooRealVar::mean = 5
0x70cf4d0/V- RooRealVar::sigma1 = 0.5
0x7817920/V- RooRealVar::sig1frac = 0.8
0x78bb070/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x73834a0/V- RooRealVar::x = 5
0x72df100/V- RooRealVar::mean = 5
0x334c6e0/V- RooRealVar::sigma2 = 1
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf206_treevistools.py.