Plot a PDF in disjunct ranges, and get normalisation right.
Usually, when comparing a fit to data, one should first plot the data, and then the PDF. In this case, the PDF is automatically normalised to match the number of data events in the plot. However, when plotting only a sub-range, when e.g. a signal region has to be blinded, one has to exclude the blinded region from the computation of the normalisation.
In this tutorial, we show how to explicitly choose the normalisation when plotting using NormRange()
.
import ROOT
x = ROOT.RooRealVar("x", "The observable", 1, 30)
tau = ROOT.RooRealVar("tau", "The exponent", -0.1337, -10.0, -0.1)
expo = ROOT.RooExponential("expo", "A falling exponential function", x, tau)
x.setRange("full", 1, 30)
x.setRange("left", 1, 10)
x.setRange("right", 20, 30)
data = expo.generate(x, 1000)
blindedData = data.reduce(CutRange="left,right")
tau.setVal(-2.0)
expo.fitTo(blindedData, Range="left,right", PrintLevel=-1)
expo.removeStringAttribute("fitrange")
canvas = ROOT.TCanvas("canvas", "canvas", 800, 600)
canvas.Divide(2, 1)
print("Now plotting with unique normalisation for each slice.\n")
canvas.cd(1)
plotFrame = x.frame(Title="Wrong: Each slice normalised over its plotting range")
blindedData.plotOn(plotFrame)
expo.plotOn(plotFrame, LineColor="r", Range="full")
expo.plotOn(plotFrame, LineColor="b", Range="left")
expo.plotOn(plotFrame, LineColor="g", Range="right")
plotFrame.Draw()
print("\n\nNow plotting with correct norm ranges:\n")
canvas.cd(2)
plotFrameWithNormRange = x.frame(Title="Right: All slices have common normalisation")
blindedData.plotOn(plotFrameWithNormRange)
expo.plotOn(plotFrameWithNormRange, LineColor="b", Range="left", NormRange="left,right")
expo.plotOn(plotFrameWithNormRange, LineColor="g", Range="right", NormRange="left,right")
expo.plotOn(plotFrameWithNormRange, LineColor="r", Range="full", NormRange="left,right", LineStyle=10)
plotFrameWithNormRange.Draw()
canvas.Draw()
canvas.SaveAs("rf212_plottingInRanges_blinding.png")
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'full' created with bounds [1,30]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'left' created with bounds [1,10]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'right' created with bounds [20,30]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_expo_expoData_left' created with bounds [1,10]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_expo_expoData_right' created with bounds [20,30]
[#1] INFO:Fitting -- RooAbsPdf::fitTo(expo_over_expo_Int[x|left,right]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_expo_over_expo_Int[x|left,right]_expoData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'full', curve is normalized to data in given range
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'left', curve is normalized to data in given range
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'right', curve is normalized to data in given range
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) only plotting range 'full'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(expo) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
Now plotting with unique normalisation for each slice.
Now plotting with correct norm ranges:
- Date
- June 2021
- Author
- Harshal Shende, Stephan Hageboeck (C++ version)
Definition in file rf212_plottingInRanges_blinding.py.