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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: complete example with use of conditional 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 fake 'external' data with per-event error
# ------------------------------------------------------------------------------------------------------
# Use landau pdf to get somewhat realistic distribution with long tail
pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25)
expDataDterr = pdfDtErr.generate({dterr}, 10000)
# Sample data from conditional decay_gm(dt|dterr)
# ---------------------------------------------------------------------------------------------
# Specify external dataset with dterr values to use decay_dm as
# conditional pdf
data = decay_gm.generate({dt}, ProtoData=expDataDterr)
# Fit conditional decay_dm(dt|dterr)
# ---------------------------------------------------------------------
# Specify dterr as conditional observable
decay_gm.fitTo(data, ConditionalObservables={dterr}, PrintLevel=-1)
# Plot conditional decay_dm(dt|dterr)
# ---------------------------------------------------------------------
# Make two-dimensional plot of conditional pdf in (dt,dterr)
hh_decay = decay_gm.createHistogram("hh_decay", dt, Binning=50, YVar=dict(var=dterr, Binning=50))
hh_decay.SetLineColor(ROOT.kBlue)
# Plot decay_gm(dt|dterr) at various values of dterr
frame = dt.frame(Title="Slices of decay(dt|dterr) at various dterr")
for ibin in range(0, 100, 20):
dterr.setBin(ibin)
decay_gm.plotOn(frame, Normalization=5.0)
# Make projection of data an dt
frame2 = dt.frame(Title="Projection of decay(dt|dterr) on dt")
data.plotOn(frame2)
# Make projection of decay(dt|dterr) on dt.
#
# Instead of integrating out dterr, a weighted average of curves
# at values dterr_i as given in the external dataset.
# (The kTRUE argument bins the data before projection to speed up the process)
decay_gm.plotOn(frame2, ProjWData=(expDataDterr, True))
# Draw all frames on canvas
c = ROOT.TCanvas("rf306_condpereventerrors", "rf306_condperventerrors", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.20)
hh_decay.GetZaxis().SetTitleOffset(2.5)
hh_decay.Draw("surf")
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.6)
frame.Draw()
c.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.6)
frame2.Draw()
c.SaveAs("rf306_condpereventerrors.png")
[#1] INFO:Fitting -- RooAbsPdf::fitTo(gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_over_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt]) 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_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_over_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt]_decay_gmData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt,dterr]) using numeric integrator RooIntegrator1D to calculate Int(dterr)
[#1] INFO:Plotting -- RooAbsReal::plotOn(decay_gm) plot on dt averages using data variables (dterr)
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf306_condpereventerrors.py.