Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf610_visualerror.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Likelihood and minimization: visualization of errors from a covariance matrix
5///
6/// \macro_image
7/// \macro_code
8/// \macro_output
9///
10/// \date April 2009
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataHist.h"
15#include "RooGaussian.h"
16#include "RooAddPdf.h"
17#include "RooPlot.h"
18#include "TCanvas.h"
19#include "TAxis.h"
20#include "TAxis.h"
21using namespace RooFit;
22
24{
25 // S e t u p e x a m p l e f i t
26 // ---------------------------------------
27
28 // Create sum of two Gaussians pdf with factory
29 RooRealVar x("x", "x", -10, 10);
30
31 RooRealVar m("m", "m", 0, -10, 10);
32 RooRealVar s("s", "s", 2, 1, 50);
33 RooGaussian sig("sig", "sig", x, m, s);
34
35 RooRealVar m2("m2", "m2", -1, -10, 10);
36 RooRealVar s2("s2", "s2", 6, 1, 50);
37 RooGaussian bkg("bkg", "bkg", x, m2, s2);
38
39 RooRealVar fsig("fsig", "fsig", 0.33, 0, 1);
40 RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
41
42 // Create binned dataset
43 x.setBins(25);
44 std::unique_ptr<RooAbsData> d{model.generateBinned(x, 1000)};
45
46 // Perform fit and save fit result
47 std::unique_ptr<RooFitResult> r{model.fitTo(*d, Save(), PrintLevel(-1))};
48
49 // V i s u a l i z e f i t e r r o r
50 // -------------------------------------
51
52 // Make plot frame
53 RooPlot *frame = x.frame(Bins(40), Title("P.d.f with visualized 1-sigma error band"));
54 d->plotOn(frame);
55
56 // Visualize 1-sigma error encoded in fit result 'r' as orange band using linear error propagation
57 // This results in an error band that is by construction symmetric
58 //
59 // The linear error is calculated as
60 // error(x) = Z* F_a(x) * Corr(a,a') F_a'(x)
61 //
62 // where F_a(x) = [ f(x,a+da) - f(x,a-da) ] / 2,
63 //
64 // with f(x) = the plotted curve
65 // 'da' = error taken from the fit result
66 // Corr(a,a') = the correlation matrix from the fit result
67 // Z = requested significance 'Z sigma band'
68 //
69 // The linear method is fast (required 2*N evaluations of the curve, where N is the number of parameters),
70 // but may not be accurate in the presence of strong correlations (~>0.9) and at Z>2 due to linear and
71 // Gaussian approximations made
72 //
73 model.plotOn(frame, VisualizeError(*r, 1), FillColor(kOrange));
74
75 // Calculate error using sampling method and visualize as dashed red line.
76 //
77 // In this method a number of curves is calculated with variations of the parameter values, as sampled
78 // from a multi-variate Gaussian pdf that is constructed from the fit results covariance matrix.
79 // The error(x) is determined by calculating a central interval that capture N% of the variations
80 // for each value of x, where N% is controlled by Z (i.e. Z=1 gives N=68%). The number of sampling curves
81 // is chosen to be such that at least 100 curves are expected to be outside the N% interval, and is minimally
82 // 100 (e.g. Z=1->Ncurve=356, Z=2->Ncurve=2156)) Intervals from the sampling method can be asymmetric,
83 // and may perform better in the presence of strong correlations, but may take (much) longer to calculate
84 model.plotOn(frame, VisualizeError(*r, 1, false), DrawOption("L"), LineWidth(2), LineColor(kRed));
85
86 // Perform the same type of error visualization on the background component only.
87 // The VisualizeError() option can generally applied to _any_ kind of plot (components, asymmetries, efficiencies
88 // etc..)
89 model.plotOn(frame, VisualizeError(*r, 1), FillColor(kOrange), Components("bkg"));
90 model.plotOn(frame, VisualizeError(*r, 1, false), DrawOption("L"), LineWidth(2), LineColor(kRed), Components("bkg"),
92
93 // Overlay central value
94 model.plotOn(frame);
95 model.plotOn(frame, Components("bkg"), LineStyle(kDashed));
96 d->plotOn(frame);
97 frame->SetMinimum(0);
98
99 // V i s u a l i z e p a r t i a l f i t e r r o r
100 // ------------------------------------------------------
101
102 // Make plot frame
103 RooPlot *frame2 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from (m,m2)"));
104
105 // Visualize partial error. For partial error visualization the covariance matrix is first reduced as follows
106 // ___ -1
107 // Vred = V22 = V11 - V12 * V22 * V21
108 //
109 // Where V11,V12,V21,V22 represent a block decomposition of the covariance matrix into observables that
110 // are propagated (labeled by index '1') and that are not propagated (labeled by index '2'), and V22bar
111 // is the Shur complement of V22, calculated as shown above
112 //
113 // (Note that Vred is _not_ a simple sub-matrix of V)
114
115 // Propagate partial error due to shape parameters (m,m2) using linear and sampling method
116 model.plotOn(frame2, VisualizeError(*r, RooArgSet(m, m2), 2), FillColor(kCyan));
117 model.plotOn(frame2, Components("bkg"), VisualizeError(*r, RooArgSet(m, m2), 2), FillColor(kCyan));
118
119 model.plotOn(frame2);
120 model.plotOn(frame2, Components("bkg"), LineStyle(kDashed));
121 frame2->SetMinimum(0);
122
123 // Make plot frame
124 RooPlot *frame3 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from (s,s2)"));
125
126 // Propagate partial error due to yield parameter using linear and sampling method
127 model.plotOn(frame3, VisualizeError(*r, RooArgSet(s, s2), 2), FillColor(kGreen));
128 model.plotOn(frame3, Components("bkg"), VisualizeError(*r, RooArgSet(s, s2), 2), FillColor(kGreen));
129
130 model.plotOn(frame3);
131 model.plotOn(frame3, Components("bkg"), LineStyle(kDashed));
132 frame3->SetMinimum(0);
133
134 // Make plot frame
135 RooPlot *frame4 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from fsig"));
136
137 // Propagate partial error due to yield parameter using linear and sampling method
138 model.plotOn(frame4, VisualizeError(*r, RooArgSet(fsig), 2), FillColor(kMagenta));
139 model.plotOn(frame4, Components("bkg"), VisualizeError(*r, RooArgSet(fsig), 2), FillColor(kMagenta));
140
141 model.plotOn(frame4);
142 model.plotOn(frame4, Components("bkg"), LineStyle(kDashed));
143 frame4->SetMinimum(0);
144
145 TCanvas *c = new TCanvas("rf610_visualerror", "rf610_visualerror", 800, 800);
146 c->Divide(2, 2);
147 c->cd(1);
148 gPad->SetLeftMargin(0.15);
149 frame->GetYaxis()->SetTitleOffset(1.4);
150 frame->Draw();
151 c->cd(2);
152 gPad->SetLeftMargin(0.15);
153 frame2->GetYaxis()->SetTitleOffset(1.6);
154 frame2->Draw();
155 c->cd(3);
156 gPad->SetLeftMargin(0.15);
157 frame3->GetYaxis()->SetTitleOffset(1.6);
158 frame3->Draw();
159 c->cd(4);
160 gPad->SetLeftMargin(0.15);
161 frame4->GetYaxis()->SetTitleOffset(1.6);
162 frame4->Draw();
163}
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
@ kOrange
Definition Rtypes.h:67
@ kGreen
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kCyan
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
#define gPad
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
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:185
virtual void SetMinimum(double minimum=-1111)
Set minimum value of Y axis.
Definition RooPlot.cxx:1007
TAxis * GetYaxis() const
Definition RooPlot.cxx:1224
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Components(Args_t &&... argsOrArgSet)
RooCmdArg DrawOption(const char *opt)
RooCmdArg FillColor(Color_t color)
RooCmdArg LineWidth(Width_t width)
RooCmdArg VisualizeError(const RooDataSet &paramData, double Z=1)
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 CodegenImpl.h:64
const char * Title
Definition TXMLSetup.cxx:68
TMarker m
Definition textangle.C:8