'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #602
Setting up a chi^2 fit to a binned dataset
from __future__ import print_function
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])
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
d = model.generate({x}, 10000)
dh = d.binnedClone()
ll = ROOT.RooLinkedList()
model.chi2FitTo(dh, ll)
dsmall = d.reduce(ROOT.RooFit.EventRange(1, 100))
dhsmall = dsmall.binnedClone()
chi2_lowstat = model.createChi2(dhsmall)
print(chi2_lowstat.getVal())
[#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: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)
Minuit2Minimizer: Minimize with max-calls 2000 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = 104.639633447510988
Edm = 0.000778057047730882148
Nfcn = 70
a0 = 0.501526 +/- 0.0229096 (limited)
a1 = 0.158456 +/- 0.0368354 (limited)
bkgfrac = 0.506609 +/- 0.011349 (limited)
sig1frac = 0.815448 +/- 0.0373695 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
90.86495991394004
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C version)
Definition in file rf602_chi2fit.py.