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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: usage of full pdf with per-event errors

import ROOT
# B-physics pdf with per-event Gaussian resolution
# ----------------------------------------------------------------------------------------------
# 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")
# Construct empirical pdf for per-event error
# -----------------------------------------------------------------
# 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 conditional product decay_dm(dt|dterr)*pdf(dterr)
# ----------------------------------------------------------------------------------------------------------------------
# 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))
# Sample, fit and plot product model
# ------------------------------------------------------------------
# Specify external dataset with dterr values to use model_dm as
# conditional pdf
data = model.generate({dt, dterr}, 10000)
# Fit conditional decay_dm(dt|dterr)
# ---------------------------------------------------------------------
# Specify dterr as conditional observable
model.fitTo(data, PrintLevel=-1)
# Plot conditional decay_dm(dt|dterr)
# ---------------------------------------------------------------------
# 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)
# 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")
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf307_fullpereventerrors.py.