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()
.
Thanks to Marc Escalier for asking how to do this correctly.
Author: Stephan Hageboeck
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, April 17, 2024 at 11:18 AM.
%%cpp -d
#include <RooDataSet.h>
#include <RooExponential.h>
#include <RooPlot.h>
#include <RooRealVar.h>
#include <TCanvas.h>
using namespace RooFit;
Make a fit model
RooRealVar x("x", "The observable", 1, 30);
RooRealVar tau("tau", "The exponent", -0.1337, -10., -0.1);
RooExponential expo("expo", "A falling exponential function", x, tau);
Define the sidebands (e.g. background regions)
x.setRange("full", 1, 30);
x.setRange("left", 1, 10);
x.setRange("right", 20, 30);
[#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]
Generate toy data, and cut out the blinded region.
std::unique_ptr<RooDataSet> data{expo.generate(x, 1000)};
std::unique_ptr<RooAbsData> blindedData{data->reduce(CutRange("left,right"))};
input_line_50:2:2: warning: 'data' shadows a declaration with the same name in the 'std' namespace; use '::data' to reference this declaration std::unique_ptr<RooDataSet> data{expo.generate(x, 1000)}; ^
The fit should be done only in the unblinded regions, otherwise it would try to make the model adapt to the empty bins in the blinded region.
tau.setVal(-2.);
expo.fitTo(*blindedData, Range("left,right"), PrintLevel(-1));
[#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
Clear the "fitrange" attribute of the PDF. Otherwise, the fitrange would be automatically taken as the NormRange() for plotting. We want to avoid this, because the point of this tutorial is to show what can go wrong when the NormRange() is not specified.
expo.removeStringAttribute("fitrange");
Here we will plot the results
TCanvas *canvas=new TCanvas("canvas","canvas",800,600);
canvas->Divide(2,1);
Plotting each slice on its own normalises the PDF over its plotting range. For the full curve, that means that the blinded region where data is missing is included in the normalisation calculation. The PDF therefore comes out too low, and doesn't match up with the slices in the side bands, which are normalised to "their" data.
std::cout << "Now plotting with unique normalisation for each slice." << std::endl;
canvas->cd(1);
RooPlot* plotFrame = x.frame(RooFit::Title("Wrong: Each slice normalised over its plotting range"));
Now plotting with unique normalisation for each slice.
Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
blindedData->plotOn(plotFrame);
expo.plotOn(plotFrame, LineColor(kRed), Range("full"));
expo.plotOn(plotFrame, LineColor(kBlue), Range("left"));
expo.plotOn(plotFrame, LineColor(kGreen), Range("right"));
plotFrame->Draw();
[#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
Make the same plot, but normalise each piece with respect to the regions "left" AND "right". This requires setting a "NormRange", which tells RooFit over which range the PDF has to be integrated to normalise. This means that the normalisation of the blue and green curves is slightly different from the left plot, because they get a common scale factor.
std::cout << "\n\nNow plotting with correct norm ranges:" << std::endl;
canvas->cd(2);
RooPlot* plotFrameWithNormRange = x.frame(RooFit::Title("Right: All slices have common normalisation"));
Now plotting with correct norm ranges:
Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
blindedData->plotOn(plotFrameWithNormRange);
expo.plotOn(plotFrameWithNormRange, LineColor(kBlue), Range("left"), RooFit::NormRange("left,right"));
expo.plotOn(plotFrameWithNormRange, LineColor(kGreen), Range("right"), RooFit::NormRange("left,right"));
expo.plotOn(plotFrameWithNormRange, LineColor(kRed), Range("full"), RooFit::NormRange("left,right"), LineStyle(10));
plotFrameWithNormRange->Draw();
canvas->Draw();
[#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'
Draw all canvases
%jsroot on
gROOT->GetListOfCanvases()->Draw()