Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf501_simultaneouspdf.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Organization and simultaneous fits: using simultaneous pdfs to describe simultaneous
5## fits to multiple datasets
6##
7## \macro_image
8## \macro_code
9## \macro_output
10##
11## \date February 2018
12## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14import ROOT
15
16
17# Create model for physics sample
18# -------------------------------------------------------------
19
20# Create observables
21x = ROOT.RooRealVar("x", "x", -8, 8)
22
23# Construct signal pdf
24mean = ROOT.RooRealVar("mean", "mean", 0, -8, 8)
25sigma = ROOT.RooRealVar("sigma", "sigma", 0.3, 0.1, 10)
26gx = ROOT.RooGaussian("gx", "gx", x, mean, sigma)
27
28# Construct background pdf
29a0 = ROOT.RooRealVar("a0", "a0", -0.1, -1, 1)
30a1 = ROOT.RooRealVar("a1", "a1", 0.004, -1, 1)
31px = ROOT.RooChebychev("px", "px", x, [a0, a1])
32
33# Construct composite pdf
34f = ROOT.RooRealVar("f", "f", 0.2, 0.0, 1.0)
35model = ROOT.RooAddPdf("model", "model", [gx, px], [f])
36
37# Create model for control sample
38# --------------------------------------------------------------
39
40# Construct signal pdf.
41# NOTE that sigma is shared with the signal sample model
42mean_ctl = ROOT.RooRealVar("mean_ctl", "mean_ctl", -3, -8, 8)
43gx_ctl = ROOT.RooGaussian("gx_ctl", "gx_ctl", x, mean_ctl, sigma)
44
45# Construct the background pdf
46a0_ctl = ROOT.RooRealVar("a0_ctl", "a0_ctl", -0.1, -1, 1)
47a1_ctl = ROOT.RooRealVar("a1_ctl", "a1_ctl", 0.5, -0.1, 1)
48px_ctl = ROOT.RooChebychev("px_ctl", "px_ctl", x, [a0_ctl, a1_ctl])
49
50# Construct the composite model
51f_ctl = ROOT.RooRealVar("f_ctl", "f_ctl", 0.5, 0.0, 1.0)
52model_ctl = ROOT.RooAddPdf("model_ctl", "model_ctl", [gx_ctl, px_ctl], [f_ctl])
53
54# Generate events for both samples
55# ---------------------------------------------------------------
56
57# Generate 1000 events in x and y from model
58data = model.generate({x}, 1000)
59data_ctl = model_ctl.generate({x}, 2000)
60
61# Create index category and join samples
62# ---------------------------------------------------------------------------
63
64# Define category to distinguish physics and control samples events
65sample = ROOT.RooCategory("sample", "sample")
66sample.defineType("physics")
67sample.defineType("control")
68
69# Construct combined dataset in (x,sample)
70combData = ROOT.RooDataSet(
71 "combData",
72 "combined data",
73 {x},
74 Index=sample,
75 Import={"physics": data, "control": data_ctl},
76)
77
78# Construct a simultaneous pdf in (x, sample)
79# -----------------------------------------------------------------------------------
80
81# Construct a simultaneous pdf using category sample as index: associate model
82# with the physics state and model_ctl with the control state
83simPdf = ROOT.RooSimultaneous("simPdf", "simultaneous pdf", {"physics": model, "control": model_ctl}, sample)
84
85# Perform a simultaneous fit
86# ---------------------------------------------------
87
88# Perform simultaneous fit of model to data and model_ctl to data_ctl
89fitResult = simPdf.fitTo(combData, PrintLevel=-1, Save=True)
90fitResult.Print()
91
92# Plot model slices on data slices
93# ----------------------------------------------------------------
94
95# Make a frame for the physics sample
96frame1 = x.frame(Title="Physics sample")
97
98# Plot all data tagged as physics sample
99combData.plotOn(frame1, Cut="sample==sample::physics")
100
101# Plot "physics" slice of simultaneous pdf.
102# NB: You *must* project the sample index category with data using ProjWData as
103# a RooSimultaneous makes no prediction on the shape in the index category and
104# can thus not be integrated. In other words: Since the PDF doesn't know the
105# number of events in the different category states, it doesn't know how much
106# of each component it has to project out. This info is read from the data.
107simPdf.plotOn(frame1, Slice=(sample, "physics"), ProjWData=(sample, combData))
108simPdf.plotOn(frame1, Slice=(sample, "physics"), Components="px", ProjWData=(sample, combData), LineStyle="--")
109
110# The same plot for the control sample slice. We do this with a different
111# approach this time, for illustration purposes. Here, we are slicing the
112# dataset and then use the data slice for the projection, because then the
113# RooFit::Slice() becomes unnecessary. This approach is more general,
114# because you can plot sums of slices by using logical or in the Cut()
115# command.
116frame2 = x.frame(Title="Control sample")
117slicedData = combData.reduce(Cut="sample==sample::control")
118slicedData.plotOn(frame2)
119simPdf.plotOn(frame2, ProjWData=(sample, slicedData))
120simPdf.plotOn(frame2, Components="px_ctl", ProjWData=(sample, slicedData), LineStyle="--")
121
122# The same plot for all the phase space. Here, we can just use the original
123# combined dataset.
124frame3 = x.frame(Title="Both samples")
125combData.plotOn(frame3)
126simPdf.plotOn(frame3, ProjWData=(sample, combData))
127simPdf.plotOn(frame3, Components="px,px_ctl", ProjWData=(sample, combData), LineStyle="--")
128
129c = ROOT.TCanvas("rf501_simultaneouspdf", "rf501_simultaneouspdf", 1200, 400)
130c.Divide(3)
131
132
133def draw(i, frame):
134 c.cd(i)
135 ROOT.gPad.SetLeftMargin(0.15)
136 frame.GetYaxis().SetTitleOffset(1.4)
137 frame.Draw()
138
139
140draw(1, frame1)
141draw(2, frame2)
142draw(3, frame3)
143
144c.SaveAs("rf501_simultaneouspdf.png")