Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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_output
8 /// \macro_code
9 /// \author 04/2009 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataHist.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooAddPdf.h"
16 #include "RooPlot.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "TAxis.h"
20 using namespace RooFit;
21 
22 void rf610_visualerror()
23 {
24  // S e t u p e x a m p l e f i t
25  // ---------------------------------------
26 
27  // Create sum of two Gaussians p.d.f. with factory
28  RooRealVar x("x", "x", -10, 10);
29 
30  RooRealVar m("m", "m", 0, -10, 10);
31  RooRealVar s("s", "s", 2, 1, 50);
32  RooGaussian sig("sig", "sig", x, m, s);
33 
34  RooRealVar m2("m2", "m2", -1, -10, 10);
35  RooRealVar s2("s2", "s2", 6, 1, 50);
36  RooGaussian bkg("bkg", "bkg", x, m2, s2);
37 
38  RooRealVar fsig("fsig", "fsig", 0.33, 0, 1);
39  RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
40 
41  // Create binned dataset
42  x.setBins(25);
43  RooAbsData *d = model.generateBinned(x, 1000);
44 
45  // Perform fit and save fit result
46  RooFitResult *r = model.fitTo(*d, Save());
47 
48  // V i s u a l i z e f i t e r r o r
49  // -------------------------------------
50 
51  // Make plot frame
52  RooPlot *frame = x.frame(Bins(40), Title("P.d.f with visualized 1-sigma error band"));
53  d->plotOn(frame);
54 
55  // Visualize 1-sigma error encoded in fit result 'r' as orange band using linear error propagation
56  // This results in an error band that is by construction symmetric
57  //
58  // The linear error is calculated as
59  // error(x) = Z* F_a(x) * Corr(a,a') F_a'(x)
60  //
61  // where F_a(x) = [ f(x,a+da) - f(x,a-da) ] / 2,
62  //
63  // with f(x) = the plotted curve
64  // 'da' = error taken from the fit result
65  // Corr(a,a') = the correlation matrix from the fit result
66  // Z = requested significance 'Z sigma band'
67  //
68  // The linear method is fast (required 2*N evaluations of the curve, where N is the number of parameters),
69  // but may not be accurate in the presence of strong correlations (~>0.9) and at Z>2 due to linear and
70  // Gaussian approximations made
71  //
72  model.plotOn(frame, VisualizeError(*r, 1), FillColor(kOrange));
73 
74  // Calculate error using sampling method and visualize as dashed red line.
75  //
76  // In this method a number of curves is calculated with variations of the parameter values, as sampled
77  // from a multi-variate Gaussian p.d.f. that is constructed from the fit results covariance matrix.
78  // The error(x) is determined by calculating a central interval that capture N% of the variations
79  // for each value of x, where N% is controlled by Z (i.e. Z=1 gives N=68%). The number of sampling curves
80  // is chosen to be such that at least 100 curves are expected to be outside the N% interval, and is minimally
81  // 100 (e.g. Z=1->Ncurve=356, Z=2->Ncurve=2156)) Intervals from the sampling method can be asymmetric,
82  // and may perform better in the presence of strong correlations, but may take (much) longer to calculate
83  model.plotOn(frame, VisualizeError(*r, 1, kFALSE), DrawOption("L"), LineWidth(2), LineColor(kRed));
84 
85  // Perform the same type of error visualization on the background component only.
86  // The VisualizeError() option can generally applied to _any_ kind of plot (components, asymmetries, efficiencies
87  // etc..)
88  model.plotOn(frame, VisualizeError(*r, 1), FillColor(kOrange), Components("bkg"));
89  model.plotOn(frame, VisualizeError(*r, 1, kFALSE), DrawOption("L"), LineWidth(2), LineColor(kRed), Components("bkg"),
90  LineStyle(kDashed));
91 
92  // Overlay central value
93  model.plotOn(frame);
94  model.plotOn(frame, Components("bkg"), LineStyle(kDashed));
95  d->plotOn(frame);
96  frame->SetMinimum(0);
97 
98  // V i s u a l i z e p a r t i a l f i t e r r o r
99  // ------------------------------------------------------
100 
101  // Make plot frame
102  RooPlot *frame2 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from (m,m2)"));
103 
104  // Visualize partial error. For partial error visualization the covariance matrix is first reduced as follows
105  // ___ -1
106  // Vred = V22 = V11 - V12 * V22 * V21
107  //
108  // Where V11,V12,V21,V22 represent a block decomposition of the covariance matrix into observables that
109  // are propagated (labeled by index '1') and that are not propagated (labeled by index '2'), and V22bar
110  // is the Shur complement of V22, calculated as shown above
111  //
112  // (Note that Vred is _not_ a simple sub-matrix of V)
113 
114  // Propagate partial error due to shape parameters (m,m2) using linear and sampling method
115  model.plotOn(frame2, VisualizeError(*r, RooArgSet(m, m2), 2), FillColor(kCyan));
116  model.plotOn(frame2, Components("bkg"), VisualizeError(*r, RooArgSet(m, m2), 2), FillColor(kCyan));
117 
118  model.plotOn(frame2);
119  model.plotOn(frame2, Components("bkg"), LineStyle(kDashed));
120  frame2->SetMinimum(0);
121 
122  // Make plot frame
123  RooPlot *frame3 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from (s,s2)"));
124 
125  // Propagate partial error due to yield parameter using linear and sampling method
126  model.plotOn(frame3, VisualizeError(*r, RooArgSet(s, s2), 2), FillColor(kGreen));
127  model.plotOn(frame3, Components("bkg"), VisualizeError(*r, RooArgSet(s, s2), 2), FillColor(kGreen));
128 
129  model.plotOn(frame3);
130  model.plotOn(frame3, Components("bkg"), LineStyle(kDashed));
131  frame3->SetMinimum(0);
132 
133  // Make plot frame
134  RooPlot *frame4 = x.frame(Bins(40), Title("Visualization of 2-sigma partial error from fsig"));
135 
136  // Propagate partial error due to yield parameter using linear and sampling method
137  model.plotOn(frame4, VisualizeError(*r, RooArgSet(fsig), 2), FillColor(kMagenta));
138  model.plotOn(frame4, Components("bkg"), VisualizeError(*r, RooArgSet(fsig), 2), FillColor(kMagenta));
139 
140  model.plotOn(frame4);
141  model.plotOn(frame4, Components("bkg"), LineStyle(kDashed));
142  frame4->SetMinimum(0);
143 
144  TCanvas *c = new TCanvas("rf610_visualerror", "rf610_visualerror", 800, 800);
145  c->Divide(2, 2);
146  c->cd(1);
147  gPad->SetLeftMargin(0.15);
148  frame->GetYaxis()->SetTitleOffset(1.4);
149  frame->Draw();
150  c->cd(2);
151  gPad->SetLeftMargin(0.15);
152  frame2->GetYaxis()->SetTitleOffset(1.6);
153  frame2->Draw();
154  c->cd(3);
155  gPad->SetLeftMargin(0.15);
156  frame3->GetYaxis()->SetTitleOffset(1.6);
157  frame3->Draw();
158  c->cd(4);
159  gPad->SetLeftMargin(0.15);
160  frame4->GetYaxis()->SetTitleOffset(1.6);
161  frame4->Draw();
162 }