Addition and convolution: tools for visualization of ROOT.RooAbsArg expression trees
Author: Clemens Lange, Wouter Verkerke (C++ version)
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, April 17, 2024 at 11:17 AM.
import ROOT
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)
[#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
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 tree to stdout
model.Print("t")
0x9779960 RooAddPdf::model = 0.602695/1 [Auto,Clean] 0x9844d70/V- RooAddPdf::bkg = 0.20539/1 [Auto,Clean] 0x980bec0/V- RooChebychev::bkg1 = 1 [Auto,Dirty] 0x50c3570/V- RooRealVar::x = 5 0x5a38480/V- RooRealVar::a0 = 0.5 0x5a2ec60/V- RooRealVar::a1 = 0 0x9879a90/V- RooRealVar::bkg1frac = 0.2 0x9831180/V- RooExponential::bkg2 = 0.00673795 [Auto,Dirty] 0x50c3570/V- RooRealVar::x = 5 0x97d5cd0/V- RooRealVar::alpha = -1 0x9860b30/V- RooRealVar::bkgfrac = 0.5 0x97840f0/V- RooAddPdf::sig = 1/1 [Auto,Clean] 0x945de70/V- RooGaussian::sig1 = 1 [Auto,Dirty] 0x50c3570/V- RooRealVar::x = 5 0x8efcc60/V- RooRealVar::mean = 5 0x8ef3570/V- RooRealVar::sigma1 = 0.5 0x840a640/V- RooRealVar::sig1frac = 0.8 0x94457a0/V- RooGaussian::sig2 = 1 [Auto,Dirty] 0x50c3570/V- RooRealVar::x = 5 0x8efcc60/V- RooRealVar::mean = 5 0x8e11130/V- RooRealVar::sigma2 = 1
Print tree to file
model.printCompactTree("", "rf206_asciitree.txt")
Print GraphViz DOT file with representation of tree
model.graphVizTree("rf206_model.dot")