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