Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf306_condpereventerrors.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: complete example with use of conditional 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 fake 'external' data with per-event error
31# ------------------------------------------------------------------------------------------------------
32
33# Use landau pdf to get somewhat realistic distribution with long tail
34pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, ROOT.RooFit.RooConst(1), ROOT.RooFit.RooConst(0.25))
35expDataDterr = pdfDtErr.generate({dterr}, 10000)
36
37# Sample data from conditional decay_gm(dt|dterr)
38# ---------------------------------------------------------------------------------------------
39
40# Specify external dataset with dterr values to use decay_dm as
41# conditional pdf
42data = decay_gm.generate({dt}, ProtoData=expDataDterr)
43
44# Fit conditional decay_dm(dt|dterr)
45# ---------------------------------------------------------------------
46
47# Specify dterr as conditional observable
48decay_gm.fitTo(data, ConditionalObservables={dterr})
49
50# Plot conditional decay_dm(dt|dterr)
51# ---------------------------------------------------------------------
52
53# Make two-dimensional plot of conditional pdf in (dt,dterr)
54hh_decay = decay_gm.createHistogram("hh_decay", dt, Binning=50, YVar=dict(var=dterr, Binning=50))
55hh_decay.SetLineColor(ROOT.kBlue)
56
57# Plot decay_gm(dt|dterr) at various values of dterr
58frame = dt.frame(Title="Slices of decay(dt|dterr) at various dterr")
59for ibin in range(0, 100, 20):
60 dterr.setBin(ibin)
61 decay_gm.plotOn(frame, Normalization=5.0)
62
63# Make projection of data an dt
64frame2 = dt.frame(Title="Projection of decay(dt|dterr) on dt")
65data.plotOn(frame2)
66
67# Make projection of decay(dt|dterr) on dt.
68#
69# Instead of integrating out dterr, a weighted average of curves
70# at values dterr_i as given in the external dataset.
71# (The kTRUE argument bins the data before projection to speed up the process)
72decay_gm.plotOn(frame2, ProjWData=(expDataDterr, True))
73
74# Draw all frames on canvas
75c = ROOT.TCanvas("rf306_condpereventerrors", "rf306_condperventerrors", 1200, 400)
76c.Divide(3)
77c.cd(1)
78ROOT.gPad.SetLeftMargin(0.20)
79hh_decay.GetZaxis().SetTitleOffset(2.5)
80hh_decay.Draw("surf")
81c.cd(2)
82ROOT.gPad.SetLeftMargin(0.15)
83frame.GetYaxis().SetTitleOffset(1.6)
84frame.Draw()
85c.cd(3)
86ROOT.gPad.SetLeftMargin(0.15)
87frame2.GetYaxis().SetTitleOffset(1.6)
88frame2.Draw()
89
90c.SaveAs("rf306_condpereventerrors.png")