Addition and convolution: setting up an extended maximum likelihood fit
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)
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])
sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
nsig = ROOT.RooRealVar("nsig", "number of signal events", 500, 0.0, 10000)
nbkg = ROOT.RooRealVar("nbkg", "number of background events", 500, 0, 10000)
model = ROOT.RooAddPdf("model", "(g1+g2)+a", [bkg, sig], [nbkg, nsig])
data = model.generate({x})
model.fitTo(data, PrintLevel=-1)
xframe = x.frame(Title="extended ML fit example")
data.plotOn(xframe)
model.plotOn(xframe, Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected))
model.plotOn(
xframe,
Components={bkg},
LineStyle=":",
Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected),
)
ras_bkg_sig2 = {bkg, sig2}
model.plotOn(
xframe,
Components=ras_bkg_sig2,
LineStyle=":",
Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected),
)
model.Print("t")
esig = ROOT.RooExtendPdf("esig", "extended signal pdf", sig, nsig)
ebkg = ROOT.RooExtendPdf("ebkg", "extended background pdf", bkg, nbkg)
model2 = ROOT.RooAddPdf("model2", "(g1+g2)+a", [ebkg, esig])
c = ROOT.TCanvas("rf202_extendedmlfit", "rf202_extendedmlfit", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.SaveAs("rf202_extendedmlfit.png")
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-1e+30, 1e+30] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-1e+30, 1e+30] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig1,sig2)
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (bkg)
[#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)
0x66f1900 RooAddPdf::model = 0.886051/1 [Auto,Clean]
0x65c10d0/V- RooChebychev::bkg = 0.733677 [Auto,Dirty]
0x5609000/V- RooRealVar::x = 5
0x61c81a0/V- RooRealVar::a0 = 0.507582 +/- 0.0795756
0x6258df0/V- RooRealVar::a1 = 0.266323 +/- 0.133724
0x667b350/V- RooRealVar::nbkg = 427.828 +/- 38.0516
0x6649cd0/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x61ee340/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x5609000/V- RooRealVar::x = 5
0x60df040/V- RooRealVar::mean = 5
0x5fb9500/V- RooRealVar::sigma1 = 0.5
0x6448920/V- RooRealVar::sig1frac = 0.641354 +/- 0.0967277
0x62570e0/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x5609000/V- RooRealVar::x = 5
0x60df040/V- RooRealVar::mean = 5
0x50e1b40/V- RooRealVar::sigma2 = 1
0x66e08d0/V- RooRealVar::nsig = 572.096 +/- 39.9008
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf202_extendedmlfit.py.