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