Logo ROOT   6.10/09
Reference Guide
rf610_visualerror.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #610
5 ///
6 /// Visualization of errors from a covariance matrix
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 04/2009 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataHist.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooAddPdf.h"
19 #include "RooPlot.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "TAxis.h"
23 using namespace RooFit ;
24 
25 
26 void rf610_visualerror()
27 {
28  // S e t u p e x a m p l e f i t
29  // ---------------------------------------
30 
31  // Create sum of two Gaussians p.d.f. with factory
32  RooRealVar x("x","x",-10,10) ;
33 
34  RooRealVar m("m","m",0,-10,10) ;
35  RooRealVar s("s","s",2,1,50) ;
36  RooGaussian sig("sig","sig",x,m,s) ;
37 
38  RooRealVar m2("m2","m2",-1,-10,10) ;
39  RooRealVar s2("s2","s2",6,1,50) ;
40  RooGaussian bkg("bkg","bkg",x,m2,s2) ;
41 
42  RooRealVar fsig("fsig","fsig",0.33,0,1) ;
43  RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
44 
45  // Create binned dataset
46  x.setBins(25) ;
47  RooAbsData* d = model.generateBinned(x,1000) ;
48 
49  // Perform fit and save fit result
50  RooFitResult* r = model.fitTo(*d,Save()) ;
51 
52 
53  // V i s u a l i z e f i t e r r o r
54  // -------------------------------------
55 
56  // Make plot frame
57  RooPlot* frame = x.frame(Bins(40),Title("P.d.f with visualized 1-sigma error band")) ;
58  d->plotOn(frame) ;
59 
60  // Visualize 1-sigma error encoded in fit result 'r' as orange band using linear error propagation
61  // This results in an error band that is by construction symmetric
62  //
63  // The linear error is calculated as
64  // error(x) = Z* F_a(x) * Corr(a,a') F_a'(x)
65  //
66  // where F_a(x) = [ f(x,a+da) - f(x,a-da) ] / 2,
67  //
68  // with f(x) = the plotted curve
69  // 'da' = error taken from the fit result
70  // Corr(a,a') = the correlation matrix from the fit result
71  // Z = requested significance 'Z sigma band'
72  //
73  // The linear method is fast (required 2*N evaluations of the curve, where N is the number of parameters),
74  // but may not be accurate in the presence of strong correlations (~>0.9) and at Z>2 due to linear and
75  // Gaussian approximations made
76  //
77  model.plotOn(frame,VisualizeError(*r,1),FillColor(kOrange)) ;
78 
79 
80  // Calculate error using sampling method and visualize as dashed red line.
81  //
82  // In this method a number of curves is calculated with variations of the parameter values, as sampled
83  // from a multi-variate Gaussian p.d.f. that is constructed from the fit results covariance matrix.
84  // The error(x) is determined by calculating a central interval that capture N% of the variations
85  // for each value of x, where N% is controlled by Z (i.e. Z=1 gives N=68%). The number of sampling curves
86  // is chosen to be such that at least 100 curves are expected to be outside the N% interval, and is minimally
87  // 100 (e.g. Z=1->Ncurve=356, Z=2->Ncurve=2156)) Intervals from the sampling method can be asymmetric,
88  // and may perform better in the presence of strong correlations, but may take (much) longer to calculate
89  model.plotOn(frame,VisualizeError(*r,1,kFALSE),DrawOption("L"),LineWidth(2),LineColor(kRed)) ;
90 
91  // Perform the same type of error visualization on the background component only.
92  // The VisualizeError() option can generally applied to _any_ kind of plot (components, asymmetries, efficiencies etc..)
93  model.plotOn(frame,VisualizeError(*r,1),FillColor(kOrange),Components("bkg")) ;
95 
96  // Overlay central value
97  model.plotOn(frame) ;
98  model.plotOn(frame,Components("bkg"),LineStyle(kDashed)) ;
99  d->plotOn(frame) ;
100  frame->SetMinimum(0) ;
101 
102 
103  // V i s u a l i z e p a r t i a l f i t e r r o r
104  // ------------------------------------------------------
105 
106  // Make plot frame
107  RooPlot* frame2 = x.frame(Bins(40),Title("Visualization of 2-sigma partial error from (m,m2)")) ;
108 
109  // Visualize partial error. For partial error visualization the covariance matrix is first reduced as follows
110  // ___ -1
111  // Vred = V22 = V11 - V12 * V22 * V21
112  //
113  // Where V11,V12,V21,V22 represent a block decomposition of the covariance matrix into observables that
114  // are propagated (labeled by index '1') and that are not propagated (labeled by index '2'), and V22bar
115  // is the Shur complement of V22, calculated as shown above
116  //
117  // (Note that Vred is _not_ a simple sub-matrix of V)
118 
119  // Propagate partial error due to shape parameters (m,m2) using linear and sampling method
120  model.plotOn(frame2,VisualizeError(*r,RooArgSet(m,m2),2),FillColor(kCyan)) ;
121  model.plotOn(frame2,Components("bkg"),VisualizeError(*r,RooArgSet(m,m2),2),FillColor(kCyan)) ;
122 
123  model.plotOn(frame2) ;
124  model.plotOn(frame2,Components("bkg"),LineStyle(kDashed)) ;
125  frame2->SetMinimum(0) ;
126 
127 
128  // Make plot frame
129  RooPlot* frame3 = x.frame(Bins(40),Title("Visualization of 2-sigma partial error from (s,s2)")) ;
130 
131  // Propagate partial error due to yield parameter using linear and sampling method
132  model.plotOn(frame3,VisualizeError(*r,RooArgSet(s,s2),2),FillColor(kGreen)) ;
133  model.plotOn(frame3,Components("bkg"),VisualizeError(*r,RooArgSet(s,s2),2),FillColor(kGreen)) ;
134 
135  model.plotOn(frame3) ;
136  model.plotOn(frame3,Components("bkg"),LineStyle(kDashed)) ;
137  frame3->SetMinimum(0) ;
138 
139 
140  // Make plot frame
141  RooPlot* frame4 = x.frame(Bins(40),Title("Visualization of 2-sigma partial error from fsig")) ;
142 
143  // Propagate partial error due to yield parameter using linear and sampling method
144  model.plotOn(frame4,VisualizeError(*r,RooArgSet(fsig),2),FillColor(kMagenta)) ;
145  model.plotOn(frame4,Components("bkg"),VisualizeError(*r,RooArgSet(fsig),2),FillColor(kMagenta)) ;
146 
147  model.plotOn(frame4) ;
148  model.plotOn(frame4,Components("bkg"),LineStyle(kDashed)) ;
149  frame4->SetMinimum(0) ;
150 
151 
152 
153  TCanvas* c = new TCanvas("rf610_visualerror","rf610_visualerror",800,800) ;
154  c->Divide(2,2) ;
155  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
156  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
157  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
158  c->cd(4) ; gPad->SetLeftMargin(0.15) ; frame4->GetYaxis()->SetTitleOffset(1.6) ; frame4->Draw() ;
159 }
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:262
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooCmdArg DrawOption(const char *opt)
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
RooCmdArg LineColor(Color_t color)
RooCmdArg VisualizeError(const RooDataSet &paramData, Double_t Z=1)
Definition: Rtypes.h:56
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:679
Definition: Rtypes.h:56
RooCmdArg Title(const char *name)
virtual void SetMinimum(Double_t minimum=-1111)
Set minimum value of Y axis.
Definition: RooPlot.cxx:959
Double_t x[n]
Definition: legend1.C:17
RooCmdArg LineStyle(Style_t style)
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
TRandom2 r(17)
TMarker * m
Definition: textangle.C:8
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:92
The Canvas class.
Definition: TCanvas.h:31
RooCmdArg FillColor(Color_t color)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg LineWidth(Width_t width)
Definition: Rtypes.h:56
RooCmdArg Save(Bool_t flag=kTRUE)
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1135
RooCmdArg Bins(Int_t nbin)
#define gPad
Definition: TVirtualPad.h:284
Definition: Rtypes.h:57
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559