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