Logo ROOT   6.10/09
Reference Guide
rf403_weightedevts.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'DATA AND CATEGORIES' RooFit tutorial macro #403
4 //
5 // Using weights in unbinned datasets
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooDataHist.h"
19 #include "RooGaussian.h"
20 #include "RooFormulaVar.h"
21 #include "RooGenericPdf.h"
22 #include "RooPolynomial.h"
23 #include "RooChi2Var.h"
24 #include "RooMinuit.h"
25 #include "TCanvas.h"
26 #include "RooPlot.h"
27 #include "RooFitResult.h"
28 using namespace RooFit ;
29 
30 
31 class TestBasic403 : public RooFitTestUnit
32 {
33 public:
34  TestBasic403(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fits with weighted datasets",refFile,writeRef,verbose) {} ;
35  Bool_t testCode() {
36 
37  // 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
38  // -------------------------------------------------------------------------------
39 
40  // Declare observable
41  RooRealVar x("x","x",-10,10) ;
42  x.setBins(40) ;
43 
44  // Construction a uniform pdf
45  RooPolynomial p0("px","px",x) ;
46 
47  // Sample 1000 events from pdf
48  RooDataSet* data = p0.generate(x,1000) ;
49 
50 
51 
52  // 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
53  // -----------------------------------------------------------------------------------
54 
55  // Construct formula to calculate (fake) weight for events
56  RooFormulaVar wFunc("w","event weight","(x*x+10)",x) ;
57 
58  // Add column with variable w to previously generated dataset
59  RooRealVar* w = (RooRealVar*) data->addColumn(wFunc) ;
60 
61  // Instruct dataset d in interpret w as event weight rather than as observable
62  data->setWeightVar(*w) ;
63 
64 
65  // 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
66  // ---------------------------------------------------------------
67 
68  // Construction quadratic polynomial pdf for fitting
69  RooRealVar a0("a0","a0",1) ;
70  RooRealVar a1("a1","a1",0,-1,1) ;
71  RooRealVar a2("a2","a2",1,0,10) ;
72  RooPolynomial p2("p2","p2",x,RooArgList(a0,a1,a2),0) ;
73 
74  // Fit quadratic polynomial to weighted data
75 
76  // NOTE: Maximum likelihood fit to weighted data does in general
77  // NOT result in correct error estimates, unless individual
78  // event weights represent Poisson statistics themselves.
79  // In general, parameter error reflect precision of SumOfWeights
80  // events rather than NumEvents events. See comparisons below
81  RooFitResult* r_ml_wgt = p2.fitTo(*data,Save()) ;
82 
83 
84 
85  // 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
86  // ---------------------------------------------------------------
87 
88  // Construct plot frame
89  RooPlot* frame = x.frame(Title("Unbinned ML fit, binned chi^2 fit to weighted data")) ;
90 
91  // Plot data using sum-of-weights-squared error rather than Poisson errors
92  data->plotOn(frame,DataError(RooAbsData::SumW2)) ;
93 
94  // Overlay result of 2nd order polynomial fit to weighted data
95  p2.plotOn(frame) ;
96 
97 
98 
99  // M L F i t o f p d f t o e q u i v a l e n t u n w e i g h t e d d a t a s e t
100  // -----------------------------------------------------------------------------------------
101 
102  // Construct a pdf with the same shape as p0 after weighting
103  RooGenericPdf genPdf("genPdf","x*x+10",x) ;
104 
105  // Sample a dataset with the same number of events as data
106  RooDataSet* data2 = genPdf.generate(x,1000) ;
107 
108  // Sample a dataset with the same number of weights as data
109  RooDataSet* data3 = genPdf.generate(x,43000) ;
110 
111  // Fit the 2nd order polynomial to both unweighted datasets and save the results for comparison
112  RooFitResult* r_ml_unw10 = p2.fitTo(*data2,Save()) ;
113  RooFitResult* r_ml_unw43 = p2.fitTo(*data3,Save()) ;
114 
115 
116  // 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
117  // ------------------------------------------------------------------------------------
118 
119  // Construct binned clone of unbinned weighted dataset
120  RooDataHist* binnedData = data->binnedClone() ;
121 
122  // Perform chi2 fit to binned weighted dataset using sum-of-weights errors
123  //
124  // NB: Within the usual approximations of a chi2 fit, a chi2 fit to weighted
125  // data using sum-of-weights-squared errors does give correct error
126  // estimates
127  RooChi2Var chi2("chi2","chi2",p2,*binnedData,DataError(RooAbsData::SumW2)) ;
128  RooMinuit m(chi2) ;
129  m.migrad() ;
130  m.hesse() ;
131 
132  // Plot chi^2 fit result on frame as well
133  RooFitResult* r_chi2_wgt = m.save() ;
134  p2.plotOn(frame,LineStyle(kDashed),LineColor(kRed),Name("p2_alt")) ;
135 
136 
137 
138  // 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
139  // ---------------------------------------------------------------------------------------------------------------
140 
141  // Note that ML fit on 1Kevt of weighted data is closer to result of ML fit on 43Kevt of unweighted data
142  // than to 1Kevt of unweighted data, whereas the reference chi^2 fit with SumW2 error gives a result closer to
143  // that of an unbinned ML fit to 1Kevt of unweighted data.
144 
145  regResult(r_ml_unw10,"rf403_ml_unw10") ;
146  regResult(r_ml_unw43,"rf403_ml_unw43") ;
147  regResult(r_ml_wgt ,"rf403_ml_wgt") ;
148  regResult(r_chi2_wgt ,"rf403_ml_chi2") ;
149  regPlot(frame,"rf403_plot1") ;
150 
151  delete binnedData ;
152  delete data ;
153  delete data2 ;
154  delete data3 ;
155 
156  return kTRUE ;
157  }
158 } ;
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)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
Double_t x[n]
Definition: legend1.C:17
RooDataHist * binnedClone(const char *newName=0, const char *newTitle=0) const
Return binned clone of this dataset.
Definition: RooDataSet.cxx:981
virtual RooAbsArg * addColumn(RooAbsArg &var, Bool_t adjustRange=kTRUE)
Add a column with the values of the given (function) argument to this dataset.
RooCmdArg DataError(Int_t)
RooCmdArg LineStyle(Style_t style)
static double p2(double t, double a, double b, double c)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooCmdArg Name(const char *name)
TMarker * m
Definition: textangle.C:8
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg Save(Bool_t flag=kTRUE)
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooMinuit is a wrapper class around TFitter/TMinuit that provides a seamless interface between the MI...
Definition: RooMinuit.h:39