Logo ROOT   6.16/01
Reference Guide
rf306_condpereventerrors.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #306
5## Complete example with use of conditional 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 fake 'external' data with per-event error
36# ------------------------------------------------------------------------------------------------------
37
38# Use landau p.d.f to get somewhat realistic 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# Sample data from conditional decay_gm(dt|dterr)
44# ---------------------------------------------------------------------------------------------
45
46# Specify external dataset with dterr values to use decay_dm as
47# conditional p.d.f.
48data = decay_gm.generate(ROOT.RooArgSet(
49 dt), ROOT.RooFit.ProtoData(expDataDterr))
50
51# Fit conditional decay_dm(dt|dterr)
52# ---------------------------------------------------------------------
53
54# Specify dterr as conditional observable
55decay_gm.fitTo(data, ROOT.RooFit.ConditionalObservables(
56 ROOT.RooArgSet(dterr)))
57
58# Plot conditional decay_dm(dt|dterr)
59# ---------------------------------------------------------------------
60
61# Make two-dimensional plot of conditional p.d.f in (dt,dterr)
62hh_decay = decay_gm.createHistogram("hh_decay", dt, ROOT.RooFit.Binning(
63 50), ROOT.RooFit.YVar(dterr, ROOT.RooFit.Binning(50)))
64hh_decay.SetLineColor(ROOT.kBlue)
65
66# Plot decay_gm(dt|dterr) at various values of dterr
67frame = dt.frame(ROOT.RooFit.Title(
68 "Slices of decay(dt|dterr) at various dterr"))
69for ibin in range(0, 100, 20):
70 dterr.setBin(ibin)
71 decay_gm.plotOn(frame, ROOT.RooFit.Normalization(5.))
72
73# Make projection of data an dt
74frame2 = dt.frame(ROOT.RooFit.Title("Projection of decay(dt|dterr) on dt"))
75data.plotOn(frame2)
76
77# Make projection of decay(dt|dterr) on dt.
78#
79# Instead of integrating out dterr, a weighted average of curves
80# at values dterr_i as given in the external dataset.
81# (The kTRUE argument bins the data before projection to speed up the process)
82decay_gm.plotOn(frame2, ROOT.RooFit.ProjWData(expDataDterr, ROOT.kTRUE))
83
84# Draw all frames on canvas
85c = ROOT.TCanvas("rf306_condpereventerrors",
86 "rf306_condperventerrors", 1200, 400)
87c.Divide(3)
88c.cd(1)
89ROOT.gPad.SetLeftMargin(0.20)
90hh_decay.GetZaxis().SetTitleOffset(2.5)
91hh_decay.Draw("surf")
92c.cd(2)
93ROOT.gPad.SetLeftMargin(0.15)
94frame.GetYaxis().SetTitleOffset(1.6)
95frame.Draw()
96c.cd(3)
97ROOT.gPad.SetLeftMargin(0.15)
98frame2.GetYaxis().SetTitleOffset(1.6)
99frame2.Draw()
100
101c.SaveAs("rf306_condpereventerrors.png")