Logo ROOT   6.16/01
Reference Guide
rf501_simultaneouspdf.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #501
6##
7## Using simultaneous p.d.f.s to describe simultaneous fits to multiple
8## datasets
9##
10## \macro_code
11##
12## \date February 2018
13## \author Clemens Lange
14## \author Wouter Verkerke (C version)
15
16
17import ROOT
18
19
20# Create model for physics sample
21# -------------------------------------------------------------
22
23# Create observables
24x = ROOT.RooRealVar("x", "x", -8, 8)
25
26# Construct signal pdf
27mean = ROOT.RooRealVar("mean", "mean", 0, -8, 8)
28sigma = ROOT.RooRealVar("sigma", "sigma", 0.3, 0.1, 10)
29gx = ROOT.RooGaussian("gx", "gx", x, mean, sigma)
30
31# Construct background pdf
32a0 = ROOT.RooRealVar("a0", "a0", -0.1, -1, 1)
33a1 = ROOT.RooRealVar("a1", "a1", 0.004, -1, 1)
34px = ROOT.RooChebychev("px", "px", x, ROOT.RooArgList(a0, a1))
35
36# Construct composite pdf
37f = ROOT.RooRealVar("f", "f", 0.2, 0., 1.)
38model = ROOT.RooAddPdf(
39 "model", "model", ROOT.RooArgList(gx, px), ROOT.RooArgList(f))
40
41# Create model for control sample
42# --------------------------------------------------------------
43
44# Construct signal pdf.
45# NOTE that sigma is shared with the signal sample model
46mean_ctl = ROOT.RooRealVar("mean_ctl", "mean_ctl", -3, -8, 8)
47gx_ctl = ROOT.RooGaussian("gx_ctl", "gx_ctl", x, mean_ctl, sigma)
48
49# Construct the background pdf
50a0_ctl = ROOT.RooRealVar("a0_ctl", "a0_ctl", -0.1, -1, 1)
51a1_ctl = ROOT.RooRealVar("a1_ctl", "a1_ctl", 0.5, -0.1, 1)
52px_ctl = ROOT.RooChebychev(
53 "px_ctl", "px_ctl", x, ROOT.RooArgList(a0_ctl, a1_ctl))
54
55# Construct the composite model
56f_ctl = ROOT.RooRealVar("f_ctl", "f_ctl", 0.5, 0., 1.)
57model_ctl = ROOT.RooAddPdf(
58 "model_ctl", "model_ctl", ROOT.RooArgList(gx_ctl, px_ctl), ROOT.RooArgList(f_ctl))
59
60# Generate events for both samples
61# ---------------------------------------------------------------
62
63# Generate 1000 events in x and y from model
64data = model.generate(ROOT.RooArgSet(x), 100)
65data_ctl = model_ctl.generate(ROOT.RooArgSet(x), 2000)
66
67# Create index category and join samples
68# ---------------------------------------------------------------------------
69
70# Define category to distinguish physics and control samples events
71sample = ROOT.RooCategory("sample", "sample")
72sample.defineType("physics")
73sample.defineType("control")
74
75# Construct combined dataset in (x,sample)
76combData = ROOT.RooDataSet("combData", "combined data", ROOT.RooArgSet(x), ROOT.RooFit.Index(
77 sample), ROOT.RooFit.Import("physics", data), ROOT.RooFit.Import("control", data_ctl))
78
79# Construct a simultaneous pdf in (x, sample)
80# -----------------------------------------------------------------------------------
81
82# Construct a simultaneous pdf using category sample as index
83simPdf = ROOT.RooSimultaneous("simPdf", "simultaneous pdf", sample)
84
85# Associate model with the physics state and model_ctl with the control
86# state
87simPdf.addPdf(model, "physics")
88simPdf.addPdf(model_ctl, "control")
89
90# Perform a simultaneous fit
91# ---------------------------------------------------
92
93# Perform simultaneous fit of model to data and model_ctl to data_ctl
94simPdf.fitTo(combData)
95
96# Plot model slices on data slices
97# ----------------------------------------------------------------
98
99# Make a frame for the physics sample
100frame1 = x.frame(ROOT.RooFit.Bins(30), ROOT.RooFit.Title("Physics sample"))
101
102# Plot all data tagged as physics sample
103combData.plotOn(frame1, ROOT.RooFit.Cut("sample==sample::physics"))
104
105# Plot "physics" slice of simultaneous pdf.
106# NBL You _must_ project the sample index category with data using ProjWData
107# as a RooSimultaneous makes no prediction on the shape in the index category
108# and can thus not be integrated
109simPdf.plotOn(frame1, ROOT.RooFit.Slice(sample, "physics"),
110 ROOT.RooFit.ProjWData(ROOT.RooArgSet(sample), combData))
111simPdf.plotOn(frame1, ROOT.RooFit.Slice(sample, "physics"), ROOT.RooFit.Components(
112 "px"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(sample), combData), ROOT.RooFit.LineStyle(ROOT.kDashed))
113
114# The same plot for the control sample slice
115frame2 = x.frame(ROOT.RooFit.Bins(30), ROOT.RooFit.Title("Control sample"))
116combData.plotOn(frame2, ROOT.RooFit.Cut("sample==sample::control"))
117simPdf.plotOn(frame2, ROOT.RooFit.Slice(sample, "control"),
118 ROOT.RooFit.ProjWData(ROOT.RooArgSet(sample), combData))
119simPdf.plotOn(frame2, ROOT.RooFit.Slice(sample, "control"), ROOT.RooFit.Components(
120 "px_ctl"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(sample), combData), ROOT.RooFit.LineStyle(ROOT.kDashed))
121
122c = ROOT.TCanvas("rf501_simultaneouspdf",
123 "rf501_simultaneouspdf", 800, 400)
124c.Divide(2)
125c.cd(1)
126ROOT.gPad.SetLeftMargin(0.15)
127frame1.GetYaxis().SetTitleOffset(1.4)
128frame1.Draw()
129c.cd(2)
130ROOT.gPad.SetLeftMargin(0.15)
131frame2.GetYaxis().SetTitleOffset(1.4)
132frame2.Draw()
133
134c.SaveAs("rf501_simultaneouspdf.png")