Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf212_plottingInRanges_blinding.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
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/// Thanks to Marc Escalier for asking how to do this correctly.
14///
15/// \macro_image
16/// \macro_code
17/// \macro_output
18///
19/// \date March 2020
20/// \author Stephan Hageboeck
21
22#include <RooDataSet.h>
23#include <RooExponential.h>
24#include <RooPlot.h>
25#include <RooRealVar.h>
26#include <TCanvas.h>
27
28using namespace RooFit;
29
31{
32 // Make a fit model
33 RooRealVar x("x", "The observable", 1, 30);
34 RooRealVar tau("tau", "The exponent", -0.1337, -10., -0.1);
35 RooExponential expo("expo", "A falling exponential function", x, tau);
36
37 // Define the sidebands (e.g. background regions)
38 x.setRange("full", 1, 30);
39 x.setRange("left", 1, 10);
40 x.setRange("right", 20, 30);
41
42 // Generate toy data, and cut out the blinded region.
43 std::unique_ptr<RooDataSet> data{expo.generate(x, 1000)};
44 std::unique_ptr<RooAbsData> blindedData{data->reduce(CutRange("left,right"))};
45
46 // Kick tau a bit, and run an unbinned fit where the blinded data are missing.
47 // ----------------------------------------------------------------------------------------------------------
48 // The fit should be done only in the unblinded regions, otherwise it would
49 // try to make the model adapt to the empty bins in the blinded region.
50 tau.setVal(-2.);
51 expo.fitTo(*blindedData, Range("left,right"), PrintLevel(-1));
52
53 // Clear the "fitrange" attribute of the PDF. Otherwise, the fitrange would
54 // be automatically taken as the NormRange() for plotting. We want to avoid
55 // this, because the point of this tutorial is to show what can go wrong when
56 // the NormRange() is not specified.
57 expo.removeStringAttribute("fitrange");
58
59
60 // Here we will plot the results
61 TCanvas *canvas=new TCanvas("canvas","canvas",800,600);
62 canvas->Divide(2,1);
63
64
65 // Wrong:
66 // ----------------------------------------------------------------------------------------------------------
67 // Plotting each slice on its own normalises the PDF over its plotting range. For the full curve, that means
68 // that the blinded region where data is missing is included in the normalisation calculation. The PDF therefore
69 // comes out too low, and doesn't match up with the slices in the side bands, which are normalised to "their" data.
70
71 std::cout << "Now plotting with unique normalisation for each slice." << std::endl;
72 canvas->cd(1);
73 RooPlot* plotFrame = x.frame(RooFit::Title("Wrong: Each slice normalised over its plotting range"));
74
75 // Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
76 blindedData->plotOn(plotFrame);
77 expo.plotOn(plotFrame, LineColor(kRed), Range("full"));
78 expo.plotOn(plotFrame, LineColor(kBlue), Range("left"));
79 expo.plotOn(plotFrame, LineColor(kGreen), Range("right"));
80
81 plotFrame->Draw();
82
83 // Right:
84 // ----------------------------------------------------------------------------------------------------------
85 // Make the same plot, but normalise each piece with respect to the regions "left" AND "right". This requires setting
86 // a "NormRange", which tells RooFit over which range the PDF has to be integrated to normalise.
87 // This means that the normalisation of the blue and green curves is slightly different from the left plot,
88 // because they get a common scale factor.
89
90 std::cout << "\n\nNow plotting with correct norm ranges:" << std::endl;
91 canvas->cd(2);
92 RooPlot* plotFrameWithNormRange = x.frame(RooFit::Title("Right: All slices have common normalisation"));
93
94 // Plot only the blinded data, and then plot the PDF over the full range as well as both sidebands
95 blindedData->plotOn(plotFrameWithNormRange);
96 expo.plotOn(plotFrameWithNormRange, LineColor(kBlue), Range("left"), RooFit::NormRange("left,right"));
97 expo.plotOn(plotFrameWithNormRange, LineColor(kGreen), Range("right"), RooFit::NormRange("left,right"));
98 expo.plotOn(plotFrameWithNormRange, LineColor(kRed), Range("full"), RooFit::NormRange("left,right"), LineStyle(10));
99
100 plotFrameWithNormRange->Draw();
101
102 canvas->Draw();
103
104}
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Exponential PDF.
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The Canvas class.
Definition TCanvas.h:23
void Draw(Option_t *option="") override
Draw a canvas.
Definition TCanvas.cxx:851
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:716
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:1153
RooCmdArg Title(const char *name)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg NormRange(const char *rangeNameList)
RooCmdArg CutRange(const char *rangeName)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
Ta Range(0, 0, 1, 1)