Logo ROOT  
Reference Guide
rf403_weightedevts.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Data and categories: using weights in unbinned datasets
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \date July 2008
11 /// \author Wouter Verkerke
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooDataHist.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooFormulaVar.h"
19 #include "RooGenericPdf.h"
20 #include "RooPolynomial.h"
21 #include "RooChi2Var.h"
22 #include "RooMinimizer.h"
23 #include "TCanvas.h"
24 #include "TAxis.h"
25 #include "RooPlot.h"
26 #include "RooFitResult.h"
27 using namespace RooFit;
28 
29 void rf403_weightedevts()
30 {
31  // C r e a t e o b s e r v a b l e a n d u n w e i g h t e d d a t a s e t
32  // -------------------------------------------------------------------------------
33 
34  // Declare observable
35  RooRealVar x("x", "x", -10, 10);
36  x.setBins(40);
37 
38  // Construction a uniform pdf
39  RooPolynomial p0("px", "px", x);
40 
41  // Sample 1000 events from pdf
42  RooDataSet *data = p0.generate(x, 1000);
43 
44  // C a l c u l a t e w e i g h t a n d m a k e d a t a s e t w e i g h t e d
45  // -----------------------------------------------------------------------------------
46 
47  // Construct formula to calculate (fake) weight for events
48  RooFormulaVar wFunc("w", "event weight", "(x*x+10)", x);
49 
50  // Add column with variable w to previously generated dataset
51  RooRealVar *w = (RooRealVar *)data->addColumn(wFunc);
52 
53  // Dataset d is now a dataset with two observable (x,w) with 1000 entries
54  data->Print();
55 
56  // Instruct dataset wdata in interpret w as event weight rather than as observable
57  RooDataSet wdata(data->GetName(), data->GetTitle(), data, *data->get(), 0, w->GetName());
58 
59  // Dataset d is now a dataset with one observable (x) with 1000 entries and a sum of weights of ~430K
60  wdata.Print();
61 
62  // U n b i n n e d M L f i t t o w e i g h t e d d a t a
63  // ---------------------------------------------------------------
64 
65  // Construction quadratic polynomial pdf for fitting
66  RooRealVar a0("a0", "a0", 1);
67  RooRealVar a1("a1", "a1", 0, -1, 1);
68  RooRealVar a2("a2", "a2", 1, 0, 10);
69  RooPolynomial p2("p2", "p2", x, RooArgList(a0, a1, a2), 0);
70 
71  // Fit quadratic polynomial to weighted data
72 
73  // NOTE: A plain Maximum likelihood fit to weighted data does in general
74  // NOT result in correct error estimates, unless individual
75  // event weights represent Poisson statistics themselves.
76  //
77  // Fit with 'wrong' errors
78  RooFitResult *r_ml_wgt = p2.fitTo(wdata, Save());
79 
80  // A first order correction to estimated parameter errors in an
81  // (unbinned) ML fit can be obtained by calculating the
82  // covariance matrix as
83  //
84  // V' = V C-1 V
85  //
86  // where V is the covariance matrix calculated from a fit
87  // to -logL = - sum [ w_i log f(x_i) ] and C is the covariance
88  // matrix calculated from -logL' = -sum [ w_i^2 log f(x_i) ]
89  // (i.e. the weights are applied squared)
90  //
91  // A fit in this mode can be performed as follows:
92 
93  RooFitResult *r_ml_wgt_corr = p2.fitTo(wdata, Save(), SumW2Error(kTRUE));
94 
95  // P l o t w e i g h e d d a t a a n d f i t r e s u l t
96  // ---------------------------------------------------------------
97 
98  // Construct plot frame
99  RooPlot *frame = x.frame(Title("Unbinned ML fit, binned chi^2 fit to weighted data"));
100 
101  // Plot data using sum-of-weights-squared error rather than Poisson errors
102  wdata.plotOn(frame, DataError(RooAbsData::SumW2));
103 
104  // Overlay result of 2nd order polynomial fit to weighted data
105  p2.plotOn(frame);
106 
107  // ML Fit of pdf to equivalent unweighted dataset
108  // -----------------------------------------------------------------------------------------
109 
110  // Construct a pdf with the same shape as p0 after weighting
111  RooGenericPdf genPdf("genPdf", "x*x+10", x);
112 
113  // Sample a dataset with the same number of events as data
114  RooDataSet *data2 = genPdf.generate(x, 1000);
115 
116  // Sample a dataset with the same number of weights as data
117  RooDataSet *data3 = genPdf.generate(x, 43000);
118 
119  // Fit the 2nd order polynomial to both unweighted datasets and save the results for comparison
120  RooFitResult *r_ml_unw10 = p2.fitTo(*data2, Save());
121  RooFitResult *r_ml_unw43 = p2.fitTo(*data3, Save());
122 
123  // C h i 2 f i t o f p d f t o b i n n e d w e i g h t e d d a t a s e t
124  // ------------------------------------------------------------------------------------
125 
126  // Construct binned clone of unbinned weighted dataset
127  RooDataHist *binnedData = wdata.binnedClone();
128  binnedData->Print("v");
129 
130  // Perform chi2 fit to binned weighted dataset using sum-of-weights errors
131  //
132  // NB: Within the usual approximations of a chi2 fit, a chi2 fit to weighted
133  // data using sum-of-weights-squared errors does give correct error
134  // estimates
135  RooChi2Var chi2("chi2", "chi2", p2, *binnedData, DataError(RooAbsData::SumW2));
136  RooMinimizer m(chi2);
137  m.migrad();
138  m.hesse();
139 
140  // Plot chi^2 fit result on frame as well
141  RooFitResult *r_chi2_wgt = m.save();
142  p2.plotOn(frame, LineStyle(kDashed), LineColor(kRed));
143 
144  // C o m p a r e f i t r e s u l t s o f c h i 2 , M L f i t s t o ( u n ) w e i g h t e d d a t a
145  // ---------------------------------------------------------------------------------------------------------------
146 
147  // Note that ML fit on 1Kevt of weighted data is closer to result of ML fit on 43Kevt of unweighted data
148  // than to 1Kevt of unweighted data, whereas the reference chi^2 fit with SumW2 error gives a result closer to
149  // that of an unbinned ML fit to 1Kevt of unweighted data.
150 
151  cout << "==> ML Fit results on 1K unweighted events" << endl;
152  r_ml_unw10->Print();
153  cout << "==> ML Fit results on 43K unweighted events" << endl;
154  r_ml_unw43->Print();
155  cout << "==> ML Fit results on 1K weighted events with a summed weight of 43K" << endl;
156  r_ml_wgt->Print();
157  cout << "==> Corrected ML Fit results on 1K weighted events with a summed weight of 43K" << endl;
158  r_ml_wgt_corr->Print();
159  cout << "==> Chi2 Fit results on 1K weighted events with a summed weight of 43K" << endl;
160  r_chi2_wgt->Print();
161 
162  new TCanvas("rf403_weightedevts", "rf403_weightedevts", 600, 600);
163  gPad->SetLeftMargin(0.15);
164  frame->GetYaxis()->SetTitleOffset(1.8);
165  frame->Draw();
166 }
RooFormulaVar.h
m
auto * m
Definition: textangle.C:8
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooAbsData::SumW2
@ SumW2
Definition: RooAbsData.h:97
RooFit::DataError
RooCmdArg DataError(Int_t)
Definition: RooGlobalFunc.cxx:157
TNamed::GetTitle
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
RooFit::SumW2Error
RooCmdArg SumW2Error(Bool_t flag)
Definition: RooGlobalFunc.cxx:210
RooArgList
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooGaussian.h
RooChi2Var
RooChi2Var implements a simple calculation from a binned dataset and a PDF.
Definition: RooChi2Var.h:25
RooAbsData::Print
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsData.h:191
x
Double_t x[n]
Definition: legend1.C:17
TCanvas.h
RooGenericPdf.h
RooDataSet.h
RooFitResult
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
RooChi2Var.h
RooDataSet::get
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
Definition: RooDataSet.cxx:1056
RooPolynomial.h
RooDataHist
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:39
RooFormulaVar
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition: RooFormulaVar.h:30
RooFit
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition: RooCFunction1Binding.h:29
RooPolynomial
RooPolynomial implements a polynomial p.d.f of the form.
Definition: RooPolynomial.h:28
RooDataHist.h
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1258
RooPlot
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
RooRealVar.h
kRed
@ kRed
Definition: Rtypes.h:66
RooFitResult.h
RooConstVar.h
RooFit::LineColor
RooCmdArg LineColor(Color_t color)
Definition: RooGlobalFunc.cxx:57
TCanvas
The Canvas class.
Definition: TCanvas.h:23
RooMinimizer.h
TAxis.h
RooGenericPdf
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
RooFitResult::Print
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:66
kDashed
@ kDashed
Definition: TAttLine.h:48
gPad
#define gPad
Definition: TVirtualPad.h:287
RooDataSet
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
TNamed::GetName
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
RooFit::LineStyle
RooCmdArg LineStyle(Style_t style)
Definition: RooGlobalFunc.cxx:58
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:37
RooFit::Save
RooCmdArg Save(Bool_t flag=kTRUE)
Definition: RooGlobalFunc.cxx:190
RooMinimizer
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Definition: RooMinimizer.h:40
RooFit::Title
RooCmdArg Title(const char *name)
Definition: RooGlobalFunc.cxx:176
RooDataSet::addColumn
virtual RooAbsArg * addColumn(RooAbsArg &var, Bool_t adjustRange=kTRUE)
Add a column with the values of the given (function) argument to this dataset.
Definition: RooDataSet.cxx:1401