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

Detailed Description

View in nbviewer Open in SWAN
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)))
import ROOT
# Setup component pdfs
# ---------------------------------------
# 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)
# 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])
# Method 1 - Two RooAddPdfs
# ------------------------------------------
# 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])
# Add signal and background
# ------------------------------------------------
# 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])
# Sample, fit and plot model
# ---------------------------------------------------
# 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)
# Plot data and PDF overlaid
xframe = x.frame(Title="Example of composite pdf=(sig1+sig2)+bkg")
data.plotOn(xframe)
model.plotOn(xframe)
# Overlay the background component of model with a dashed line
model.plotOn(xframe, Components={bkg}, LineStyle="--")
# Overlay the background+sig2 components of model with a dotted line
model.plotOn(xframe, Components={bkg, sig2}, LineStyle=":")
# Print structure of composite pdf
model.Print("t")
# Method 2 - One RooAddPdf with recursive fractions
# ---------------------------------------------------
# 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 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
# Plot recursive addition model
# ---------------------------------------------------------
model2.plotOn(xframe, LineColor="r", LineStyle="--")
model2.plotOn(xframe, Components={bkg, sig2}, LineColor="r", LineStyle="--")
model2.Print("t")
# 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")
[#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.
[#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 -mavx2
[#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
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg,sig2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: (sig)
0x96eb0e0 RooAddPdf::model = 0.886326/1 [Auto,Clean]
0x95c6050/V- RooChebychev::bkg = 0.734412 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x7748e20/V- RooRealVar::a0 = 0.506755 +/- 0.0795919
0x6a40d60/V- RooRealVar::a1 = 0.265588 +/- 0.133931
0x967f580/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013
0x91b25d0/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x9145280/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x8a454d0/V- RooRealVar::mean = 5
0x888d940/V- RooRealVar::sigma1 = 0.5
0x6a37540/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095
0x9127350/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x8a454d0/V- RooRealVar::mean = 5
0x7ce0bd0/V- RooRealVar::sigma2 = 1
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg,sig2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
0xa0e73d0 RooAddPdf::model = 0.886326/1 [Auto,Clean]
0x95c6050/V- RooChebychev::bkg = 0.734412 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x7748e20/V- RooRealVar::a0 = 0.506755 +/- 0.0795919
0x6a40d60/V- RooRealVar::a1 = 0.265588 +/- 0.133931
0x967f580/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013
0x9145280/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x8a454d0/V- RooRealVar::mean = 5
0x888d940/V- RooRealVar::sigma1 = 0.5
0xa06c8a0/V- RooRecursiveFraction::model_recursive_fraction_sig1_2 = 0.367214 [Auto,Clean]
0x6a37540/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095
0x967f580/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013
0x9127350/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x8a82800/V- RooRealVar::x = 5
0x8a454d0/V- RooRealVar::mean = 5
0x7ce0bd0/V- RooRealVar::sigma2 = 1
0xa0e2b30/V- RooRecursiveFraction::model_recursive_fraction_sig2_3 = 0.204778 [Auto,Clean]
0x7cec0c0/V- RooConstVar::1 = 1
0x6a37540/V- RooRealVar::sig1frac = 0.641992 +/- 0.0969095
0x967f580/V- RooRealVar::bkgfrac = 0.428008 +/- 0.0356013
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf201_composite.py.