Multidimensional models: usage of full pdf with per-event errors
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 Tuesday, March 19, 2024 at 07:15 PM.
import ROOT
Observables
dt = ROOT.RooRealVar("dt", "dt", -10, 10)
dterr = ROOT.RooRealVar("dterr", "per-event error on dt", 0.01, 10)
Build a gaussian resolution model scaled by the per-error = gauss(dt,bias,sigma*dterr)
bias = ROOT.RooRealVar("bias", "bias", 0, -10, 10)
sigma = ROOT.RooRealVar("sigma", "per-event error scale factor", 1, 0.1, 10)
gm = ROOT.RooGaussModel("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr)
Construct decay(dt) (x) gauss1(dt|dterr)
tau = ROOT.RooRealVar("tau", "tau", 1.548)
decay_gm = ROOT.RooDecay("decay_gm", "decay", dt, tau, gm, type="DoubleSided")
Use landau pdf to get empirical distribution with long tail
pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25)
expDataDterr = pdfDtErr.generate({dterr}, 10000)
Construct a histogram pdf to describe the shape of the dtErr distribution
expHistDterr = expDataDterr.binnedClone()
pdfErr = ROOT.RooHistPdf("pdfErr", "pdfErr", {dterr}, expHistDterr)
Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
model = ROOT.RooProdPdf("model", "model", {pdfErr}, Conditional=({decay_gm}, {dt}))
(Alternatively you could also use the landau shape pdfDtErr) ROOT.RooProdPdf model("model", "model",pdfDtErr, ROOT.RooFit.Conditional(decay_gm,dt))
Specify external dataset with dterr values to use model_dm as conditional pdf
data = model.generate({dt, dterr}, 10000)
Specify dterr as conditional observable
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 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
Make two-dimensional plot of conditional pdf in (dt,dterr)
hh_model = model.createHistogram("hh_model", dt, Binning=50, YVar=dict(var=dterr, Binning=50))
hh_model.SetLineColor(ROOT.kBlue)
Make projection of data an dt
frame = dt.frame(Title="Projection of model(dt|dterr) on dt")
data.plotOn(frame)
model.plotOn(frame)
<cppyy.gbl.RooPlot object at 0xa3de960>
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on dt integrates over variables (dterr) [#1] INFO:NumericIntegration -- RooRealIntegral::init([pdfErr_NORM[dterr]_X_decay_gm_NORM[dt]]_Int[dterr]) using numeric integrator RooIntegrator1D to calculate Int(dterr)
Draw all frames on canvas
c = ROOT.TCanvas("rf307_fullpereventerrors", "rf307_fullpereventerrors", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.20)
hh_model.GetZaxis().SetTitleOffset(2.5)
hh_model.Draw("surf")
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.6)
frame.Draw()
c.SaveAs("rf307_fullpereventerrors.png")
Info in <TCanvas::Print>: png file rf307_fullpereventerrors.png has been created
Draw all canvases
from ROOT import gROOT
gROOT.GetListOfCanvases().Draw()