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