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)
23exp = ROOT.RooExponential("exp", "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 = exp.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# ----------------------------------------------------------------------------------------------------------
36tau.setVal(-2.0)
37exp.fitTo(blindedData)
38
39
40# Here we will plot the results
41canvas = ROOT.TCanvas("canvas", "canvas", 800, 600)
42canvas.Divide(2, 1)
43
44
45# Wrong:
46# Plotting each slice on its own normalises the PDF over its plotting range. For the full curve, that means
47# that the blinded region where data is missing is included in the normalisation calculation. The PDF therefore
48# comes out too low, and doesn't match up with the slices in the side bands, which are normalised to "their" data.
49# ----------------------------------------------------------------------------------------------------------
50
51print("Now plotting with unique normalisation for each slice.\n")
52canvas.cd(1)
53plotFrame = x.frame(Title="Wrong: Each slice normalised over its plotting range")
54
55# Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
56blindedData.plotOn(plotFrame)
57exp.plotOn(plotFrame, LineColor="r", Range="full")
58exp.plotOn(plotFrame, LineColor="b", Range="left")
59exp.plotOn(plotFrame, LineColor="g", Range="right")
60
61plotFrame.Draw()
62
63# Right:
64# Make the same plot, but normalise each piece with respect to the regions "left" AND "right". This requires setting
65# a "NormRange", which tells RooFit over which range the PDF has to be integrated to normalise.
66# This is means that the normalisation of the blue and green curves is slightly different from the left plot,
67# because they get a common scale factor.
68# ----------------------------------------------------------------------------------------------------------
69
70print("\n\nNow plotting with correct norm ranges:\n")
71canvas.cd(2)
72plotFrameWithNormRange = x.frame(Title="Right: All slices have common normalisation")
73
74# Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
75blindedData.plotOn(plotFrameWithNormRange)
76exp.plotOn(plotFrameWithNormRange, LineColor="b", Range="left", NormRange="left,right")
77exp.plotOn(plotFrameWithNormRange, LineColor="g", Range="right", NormRange="left,right")
78exp.plotOn(plotFrameWithNormRange, LineColor="r", Range="full", NormRange="left,right", LineStyle=10)
79
80plotFrameWithNormRange.Draw()
81
82canvas.Draw()
83
84canvas.SaveAs("rf212_plottingInRanges_blinding.png")