Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf212_plottingInRanges_blinding.C File Reference

Detailed Description

View in nbviewer Open in SWAN 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.

#include <RooDataSet.h>
#include <RooExponential.h>
#include <RooPlot.h>
#include <RooRealVar.h>
#include <TCanvas.h>
using namespace RooFit;
void rf212_plottingInRanges_blinding()
{
// Make a fit model
RooRealVar x("x", "The observable", 1, 30);
RooRealVar tau("tau", "The exponent", -0.1337, -10., -0.1);
RooExponential exp("exp", "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);
// Generate toy data, and cut out the blinded region.
RooDataSet* data = exp.generate(x, 1000);
auto blindedData = data->reduce(CutRange("left,right"));
// Kick tau a bit, and run an unbinned fit where the blinded data are missing.
// ----------------------------------------------------------------------------------------------------------
tau.setVal(-2.);
exp.fitTo(*blindedData);
// Here we will plot the results
TCanvas *canvas=new TCanvas("canvas","canvas",800,600);
canvas->Divide(2,1);
// Wrong:
// 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"));
// Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
blindedData->plotOn(plotFrame);
exp.plotOn(plotFrame, LineColor(kRed), Range("full"));
exp.plotOn(plotFrame, LineColor(kBlue), Range("left"));
exp.plotOn(plotFrame, LineColor(kGreen), Range("right"));
plotFrame->Draw();
// Right:
// 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 is 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"));
// Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
blindedData->plotOn(plotFrameWithNormRange);
exp.plotOn(plotFrameWithNormRange, LineColor(kBlue), Range("left"), RooFit::NormRange("left,right"));
exp.plotOn(plotFrameWithNormRange, LineColor(kGreen), Range("right"), RooFit::NormRange("left,right"));
exp.plotOn(plotFrameWithNormRange, LineColor(kRed), Range("full"), RooFit::NormRange("left,right"), LineStyle(10));
plotFrameWithNormRange->Draw();
canvas->Draw();
}
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
double exp(double)
RooAbsData * reduce(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg())
Create a reduced copy of this dataset.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
Exponential PDF.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:249
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
The Canvas class.
Definition TCanvas.h:23
void Draw(Option_t *option="") override
Draw a canvas.
Definition TCanvas.cxx:845
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:708
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1177
RooCmdArg Title(const char *name)
RooCmdArg NormRange(const char *rangeNameList)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Ta Range(0, 0, 1, 1)
␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#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:Minization -- RooMinimizer::optimizeConst: activating const optimization
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 tau -2.00000e+00 9.50000e-01 -1.00000e+01 -1.00000e-01
**********
** 3 **SET ERR 0.5
**********
**********
** 4 **SET PRINT 1
**********
**********
** 5 **SET STR 1
**********
NOW USING STRATEGY 1: TRY TO BALANCE SPEED AGAINST RELIABILITY
**********
** 6 **MIGRAD 500 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=6752.91 FROM MIGRAD STATUS=INITIATE 4 CALLS 5 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 tau -2.00000e+00 9.50000e-01 2.51381e-01 -1.27116e+04
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=1941.97 FROM MIGRAD STATUS=CONVERGED 30 CALLS 31 TOTAL
EDM=3.68147e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 tau -2.04501e-01 7.81359e-03 2.33650e-04 -7.85653e-02
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 1 ERR DEF=0.5
6.105e-05
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 500
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=1941.97 FROM HESSE STATUS=OK 5 CALLS 36 TOTAL
EDM=3.6779e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 tau -2.04501e-01 7.81358e-03 4.67300e-05 1.36495e+00
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 1 ERR DEF=0.5
6.105e-05
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
Now plotting with unique normalisation for each slice.
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'full', curve is normalized to data in given range
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'left', curve is normalized to data in given range
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'right', curve is normalized to data in given range
Now plotting with correct norm ranges:
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) only plotting range 'full'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(exp) p.d.f. curve is normalized using explicit choice of ranges 'left,right'
Date
March 2020
Author
Stephan Hageboeck

Definition in file rf212_plottingInRanges_blinding.C.