Addition and convolution: composite pdf with signal and background component
pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
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 Saturday, July 12, 2025 at 03:44 PM.
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.
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)
bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
Add signal components
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])
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])
Generate a data sample of 1000 events in x from model
data = model.generate({x}, 1000)
Fit model to data
model.fitTo(data, PrintLevel=-1)
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data [#1] INFO:Fitting -- using generic CPU library compiled with no vectorizations [#1] INFO:Fitting -- Creation of NLL object took 1.22116 ms [#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
Plot data and PDF overlaid
xframe = x.frame(Title="Example of composite pdf=(sig1+sig2)+bkg")
data.plotOn(xframe)
model.plotOn(xframe)
<cppyy.gbl.RooPlot object at 0x55999e215bd0>
Overlay the background component of model with a dashed line
model.plotOn(xframe, Components={bkg}, LineStyle="--")
<cppyy.gbl.RooPlot object at 0x55999e215bd0>
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg) [#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
Overlay the background+sig2 components of model with a dotted line
model.plotOn(xframe, Components={bkg, sig2}, LineStyle=":")
<cppyy.gbl.RooPlot object at 0x55999e215bd0>
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg,sig2) [#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: (sig)
Print structure of composite pdf
model.Print("t")
0x55999dbd1200 RooAddPdf::model = 0.886326/1 [Auto,Clean] 0x55999dbeeb00/V- RooChebychev::bkg = 0.734412 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x55999d7910d0/V- RooRealVar::a0 = 0.506755 +/- 0.0795919 0x55999ca0b560/V- RooRealVar::a1 = 0.265588 +/- 0.133931 0x55999da16030/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013 0x55999d717ce0/V- RooAddPdf::sig = 1/1 [Auto,Clean] 0x55999d8ab4e0/V- RooGaussian::sig1 = 1 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x559999056c70/V- RooRealVar::mean = 5 0x55999d3097e0/V- RooRealVar::sigma1 = 0.5 0x55999d9f0290/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095 0x55999d887230/V- RooGaussian::sig2 = 1 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x559999056c70/V- RooRealVar::mean = 5 0x559997cb07f0/V- RooRealVar::sigma2 = 1
Construct sum of models on one go using recursive fraction interpretations
model2 = bkg + (sig1 + sig2)
model2 = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig1, sig2], [bkgfrac, sig1frac], True)
NB: Each coefficient is interpreted as the fraction of the left-hand component of the i-th recursive sum, i.e.
sum4 = A + ( B + ( C + D) with fraction fA, and fC expands to
sum4 = fAA + (1-fA)(fBB + (1-fB)(fCC + (1-fC)D))
model2.plotOn(xframe, LineColor="r", LineStyle="--")
model2.plotOn(xframe, Components={bkg, sig2}, LineColor="r", LineStyle="--")
model2.Print("t")
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg,sig2) [#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: () 0x55999e26ab80 RooAddPdf::model = 0.886326/1 [Auto,Clean] 0x55999dbeeb00/V- RooChebychev::bkg = 0.734412 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x55999d7910d0/V- RooRealVar::a0 = 0.506755 +/- 0.0795919 0x55999ca0b560/V- RooRealVar::a1 = 0.265588 +/- 0.133931 0x55999da16030/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013 0x55999d8ab4e0/V- RooGaussian::sig1 = 1 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x559999056c70/V- RooRealVar::mean = 5 0x55999d3097e0/V- RooRealVar::sigma1 = 0.5 0x55999e5b3f00/V- RooRecursiveFraction::model_recursive_fraction_sig1_2 = 0.367214 [Auto,Clean] 0x55999d9f0290/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095 0x55999da16030/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013 0x55999d887230/V- RooGaussian::sig2 = 1 [Auto,Dirty] 0x55999d5472b0/V- RooRealVar::x = 5 0x559999056c70/V- RooRealVar::mean = 5 0x559997cb07f0/V- RooRealVar::sigma2 = 1 0x55999e30f820/V- RooRecursiveFraction::model_recursive_fraction_sig2_3 = 0.204778 [Auto,Clean] 0x55999e4ef7e0/V- RooConstVar::1 = 1 0x55999d9f0290/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095 0x55999da16030/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013
Draw the frame on the canvas
c = ROOT.TCanvas("rf201_composite", "rf201_composite", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.SaveAs("rf201_composite.png")
Info in <TCanvas::Print>: png file rf201_composite.png has been created
Draw all canvases
from ROOT import gROOT
gROOT.GetListOfCanvases().Draw()