Multidimensional models: working with parameterized ranges in a fit.
This an example of a fit with an acceptance that changes per-event
pdf = exp(-t/tau)
with t[tmin,5]
where t
and tmin
are both observables in the dataset
import ROOT
t = ROOT.RooRealVar("t", "t", 0, 5)
tmin = ROOT.RooRealVar("tmin", "tmin", 0, 0, 5)
t.setRange(tmin, ROOT.RooFit.RooConst(t.getMax()))
tau = ROOT.RooRealVar("tau", "tau", -1.54, -10, -0.1)
model = ROOT.RooExponential("model", "model", t, tau)
dall = model.generate({t}, 10000)
tmp = ROOT.RooGaussian("gmin", "gmin", tmin, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(0.5)).generate({tmin}, 5000)
dacc = model.generate({t}, ProtoData=tmp)
r = model.fitTo(dacc, Save=True, PrintLevel=-1)
frame = t.frame(Title="Fit to data with per-event acceptance")
dall.plotOn(frame, MarkerColor="r", LineColor="r")
model.plotOn(frame)
dacc.plotOn(frame)
r.Print("v")
c = ROOT.TCanvas("rf314_paramranges", "rf314_paramranges", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.6)
frame.Draw()
c.SaveAs("rf314_paramranges.png")
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 5000 will supercede previous event count of 10000 for normalization of PDF projections
RooFitResult: minimized FCN value: 2823.97, estimated distance to minimum: 6.74739e-08
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter InitialValue FinalValue +/- Error GblCorr.
-------------------- ------------ -------------------------- --------
tau -1.5400e+00 -1.5335e+00 +/- 2.22e-02 <none>
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf314_paramfitrange.py.