Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf212_plottingInRanges_blinding.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Plot a PDF in disjunct ranges, and get normalisation right.
5##
6## Usually, when comparing a fit to data, one should first plot the data, and then the PDF.
7## In this case, the PDF is automatically normalised to match the number of data events in the plot.
8## However, when plotting only a sub-range, when e.g. a signal region has to be blinded,
9## one has to exclude the blinded region from the computation of the normalisation.
10##
11## In this tutorial, we show how to explicitly choose the normalisation when plotting using `NormRange()`.
12##
13## \macro_code
14##
15## \date June 2021
16## \author Harshal Shende, Stephan Hageboeck (C++ version)
17
18import ROOT
19
20# Make a fit model
21x = ROOT.RooRealVar("x", "The observable", 1, 30)
22tau = ROOT.RooRealVar("tau", "The exponent", -0.1337, -10.0, -0.1)
23expo = ROOT.RooExponential("expo", "A falling exponential function", x, tau)
24
25# Define the sidebands (e.g. background regions)
26x.setRange("full", 1, 30)
27x.setRange("left", 1, 10)
28x.setRange("right", 20, 30)
29
30# Generate toy data, and cut out the blinded region.
31data = expo.generate(x, 1000)
32blindedData = data.reduce(CutRange="left,right")
33
34# Kick tau a bit, and run an unbinned fit where the blinded data are missing.
35# ----------------------------------------------------------------------------------------------------------
36# The fit should be done only in the unblinded regions, otherwise it would try
37# to make the model adapt to the empty bins in the blinded region.
38tau.setVal(-2.0)
39expo.fitTo(blindedData, Range="left,right")
40
41# Clear the "fitrange" attribute of the PDF. Otherwise, the fitrange would be
42# automatically taken as the NormRange() for plotting. We want to avoid this,
43# because the point of this tutorial is to show what can go wrong when the
44# NormRange() is not specified.
45expo.removeStringAttribute("fitrange")
46
47
48# Here we will plot the results
49canvas = ROOT.TCanvas("canvas", "canvas", 800, 600)
50canvas.Divide(2, 1)
51
52
53# Wrong:
54# ----------------------------------------------------------------------------------------------------------
55# Plotting each slice on its own normalises the PDF over its plotting range. For the full curve, that means
56# that the blinded region where data is missing is included in the normalisation calculation. The PDF therefore
57# comes out too low, and doesn't match up with the slices in the side bands, which are normalised to "their" data.
58
59print("Now plotting with unique normalisation for each slice.\n")
60canvas.cd(1)
61plotFrame = x.frame(Title="Wrong: Each slice normalised over its plotting range")
62
63# Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
64blindedData.plotOn(plotFrame)
65expo.plotOn(plotFrame, LineColor="r", Range="full")
66expo.plotOn(plotFrame, LineColor="b", Range="left")
67expo.plotOn(plotFrame, LineColor="g", Range="right")
68
69plotFrame.Draw()
70
71# Right:
72# ----------------------------------------------------------------------------------------------------------
73# Make the same plot, but normalise each piece with respect to the regions "left" AND "right". This requires setting
74# a "NormRange", which tells RooFit over which range the PDF has to be integrated to normalise.
75# This means that the normalisation of the blue and green curves is slightly different from the left plot,
76# because they get a common scale factor.
77
78print("\n\nNow plotting with correct norm ranges:\n")
79canvas.cd(2)
80plotFrameWithNormRange = x.frame(Title="Right: All slices have common normalisation")
81
82# Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
83blindedData.plotOn(plotFrameWithNormRange)
84expo.plotOn(plotFrameWithNormRange, LineColor="b", Range="left", NormRange="left,right")
85expo.plotOn(plotFrameWithNormRange, LineColor="g", Range="right", NormRange="left,right")
86expo.plotOn(plotFrameWithNormRange, LineColor="r", Range="full", NormRange="left,right", LineStyle=10)
87
88plotFrameWithNormRange.Draw()
89
90canvas.Draw()
91
92canvas.SaveAs("rf212_plottingInRanges_blinding.png")