Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf307_fullpereventerrors.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: usage of full pdf with per-event errors
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# B-physics pdf with per-event Gaussian resolution
14# ----------------------------------------------------------------------------------------------
15
16# Observables
17dt = ROOT.RooRealVar("dt", "dt", -10, 10)
18dterr = ROOT.RooRealVar("dterr", "per-event error on dt", 0.01, 10)
19
20# Build a gaussian resolution model scaled by the per-error =
21# gauss(dt,bias,sigma*dterr)
22bias = ROOT.RooRealVar("bias", "bias", 0, -10, 10)
23sigma = ROOT.RooRealVar("sigma", "per-event error scale factor", 1, 0.1, 10)
24gm = ROOT.RooGaussModel("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr)
25
26# Construct decay(dt) (x) gauss1(dt|dterr)
27tau = ROOT.RooRealVar("tau", "tau", 1.548)
28decay_gm = ROOT.RooDecay("decay_gm", "decay", dt, tau, gm, type="DoubleSided")
29
30# Construct empirical pdf for per-event error
31# -----------------------------------------------------------------
32
33# Use landau pdf to get empirical distribution with long tail
34pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25)
35expDataDterr = pdfDtErr.generate({dterr}, 10000)
36
37# Construct a histogram pdf to describe the shape of the dtErr distribution
38expHistDterr = expDataDterr.binnedClone()
39pdfErr = ROOT.RooHistPdf("pdfErr", "pdfErr", {dterr}, expHistDterr)
40
41# Construct conditional product decay_dm(dt|dterr)*pdf(dterr)
42# ----------------------------------------------------------------------------------------------------------------------
43
44# Construct production of conditional decay_dm(dt|dterr) with empirical
45# pdfErr(dterr)
46model = ROOT.RooProdPdf("model", "model", {pdfErr}, Conditional=({decay_gm}, {dt}))
47
48# (Alternatively you could also use the landau shape pdfDtErr)
49# ROOT.RooProdPdf model("model", "model",pdfDtErr,
50# ROOT.RooFit.Conditional(decay_gm,dt))
51
52# Sample, fit and plot product model
53# ------------------------------------------------------------------
54
55# Specify external dataset with dterr values to use model_dm as
56# conditional pdf
57data = model.generate({dt, dterr}, 10000)
58
59# Fit conditional decay_dm(dt|dterr)
60# ---------------------------------------------------------------------
61
62# Specify dterr as conditional observable
63model.fitTo(data, PrintLevel=-1)
64
65# Plot conditional decay_dm(dt|dterr)
66# ---------------------------------------------------------------------
67
68# Make two-dimensional plot of conditional pdf in (dt,dterr)
69hh_model = model.createHistogram("hh_model", dt, Binning=50, YVar=dict(var=dterr, Binning=50))
70hh_model.SetLineColor(ROOT.kBlue)
71
72# Make projection of data an dt
73frame = dt.frame(Title="Projection of model(dt|dterr) on dt")
74data.plotOn(frame)
75model.plotOn(frame)
76
77# Draw all frames on canvas
78c = ROOT.TCanvas("rf307_fullpereventerrors", "rf307_fullpereventerrors", 800, 400)
79c.Divide(2)
80c.cd(1)
81ROOT.gPad.SetLeftMargin(0.20)
82hh_model.GetZaxis().SetTitleOffset(2.5)
83hh_model.Draw("surf")
84c.cd(2)
85ROOT.gPad.SetLeftMargin(0.15)
86frame.GetYaxis().SetTitleOffset(1.6)
87frame.Draw()
88
89c.SaveAs("rf307_fullpereventerrors.png")