rf403_weightedevts.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'DATA AND CATEGORIES' RooFit tutorial macro #403
00004 // 
00005 // Using weights in unbinned datasets
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooDataHist.h"
00019 #include "RooGaussian.h"
00020 #include "RooFormulaVar.h"
00021 #include "RooGenericPdf.h"
00022 #include "RooPolynomial.h"
00023 #include "RooChi2Var.h"
00024 #include "RooMinuit.h"
00025 #include "TCanvas.h"
00026 #include "RooPlot.h"
00027 #include "RooFitResult.h"
00028 using namespace RooFit ;
00029 
00030 
00031 class TestBasic403 : public RooFitTestUnit
00032 {
00033 public: 
00034   TestBasic403(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fits with weighted datasets",refFile,writeRef,verbose) {} ;
00035   Bool_t testCode() {
00036 
00037   // 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 
00038   // -------------------------------------------------------------------------------
00039 
00040   // Declare observable
00041   RooRealVar x("x","x",-10,10) ;
00042   x.setBins(40) ;
00043 
00044   // Construction a uniform pdf
00045   RooPolynomial p0("px","px",x) ;
00046 
00047   // Sample 1000 events from pdf
00048   RooDataSet* data = p0.generate(x,1000) ;
00049 
00050  
00051 
00052   // 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 
00053   // -----------------------------------------------------------------------------------
00054 
00055   // Construct formula to calculate (fake) weight for events
00056   RooFormulaVar wFunc("w","event weight","(x*x+10)",x) ;
00057 
00058   // Add column with variable w to previously generated dataset
00059   RooRealVar* w = (RooRealVar*) data->addColumn(wFunc) ;
00060 
00061   // Instruct dataset d in interpret w as event weight rather than as observable
00062   data->setWeightVar(*w) ;
00063 
00064 
00065   // 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 
00066   // ---------------------------------------------------------------
00067 
00068   // Construction quadratic polynomial pdf for fitting
00069   RooRealVar a0("a0","a0",1) ;
00070   RooRealVar a1("a1","a1",0,-1,1) ;
00071   RooRealVar a2("a2","a2",1,0,10) ;
00072   RooPolynomial p2("p2","p2",x,RooArgList(a0,a1,a2),0) ;
00073 
00074   // Fit quadratic polynomial to weighted data
00075 
00076   // NOTE: Maximum likelihood fit to weighted data does in general 
00077   //       NOT result in correct error estimates, unless individual
00078   //       event weights represent Poisson statistics themselves.
00079   //       In general, parameter error reflect precision of SumOfWeights
00080   //       events rather than NumEvents events. See comparisons below
00081   RooFitResult* r_ml_wgt = p2.fitTo(*data,Save()) ;
00082 
00083 
00084 
00085   // 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 
00086   // ---------------------------------------------------------------
00087 
00088   // Construct plot frame
00089   RooPlot* frame = x.frame(Title("Unbinned ML fit, binned chi^2 fit to weighted data")) ;
00090 
00091   // Plot data using sum-of-weights-squared error rather than Poisson errors
00092   data->plotOn(frame,DataError(RooAbsData::SumW2)) ;
00093 
00094   // Overlay result of 2nd order polynomial fit to weighted data
00095   p2.plotOn(frame) ;
00096 
00097 
00098 
00099   // 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
00100   // -----------------------------------------------------------------------------------------
00101   
00102   // Construct a pdf with the same shape as p0 after weighting
00103   RooGenericPdf genPdf("genPdf","x*x+10",x) ;
00104 
00105   // Sample a dataset with the same number of events as data
00106   RooDataSet* data2 = genPdf.generate(x,1000) ;
00107 
00108   // Sample a dataset with the same number of weights as data
00109   RooDataSet* data3 = genPdf.generate(x,43000) ;
00110 
00111   // Fit the 2nd order polynomial to both unweighted datasets and save the results for comparison
00112   RooFitResult* r_ml_unw10 = p2.fitTo(*data2,Save()) ;
00113   RooFitResult* r_ml_unw43 = p2.fitTo(*data3,Save()) ;
00114 
00115 
00116   // 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
00117   // ------------------------------------------------------------------------------------
00118 
00119   // Construct binned clone of unbinned weighted dataset
00120   RooDataHist* binnedData = data->binnedClone() ;
00121 
00122   // Perform chi2 fit to binned weighted dataset using sum-of-weights errors
00123   // 
00124   // NB: Within the usual approximations of a chi2 fit, a chi2 fit to weighted
00125   // data using sum-of-weights-squared errors does give correct error
00126   // estimates
00127   RooChi2Var chi2("chi2","chi2",p2,*binnedData,DataError(RooAbsData::SumW2)) ;
00128   RooMinuit m(chi2) ;
00129   m.migrad() ;
00130   m.hesse() ;
00131 
00132   // Plot chi^2 fit result on frame as well
00133   RooFitResult* r_chi2_wgt = m.save() ;
00134   p2.plotOn(frame,LineStyle(kDashed),LineColor(kRed),Name("p2_alt")) ;
00135 
00136 
00137 
00138   // 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 
00139   // ---------------------------------------------------------------------------------------------------------------
00140 
00141   // Note that ML fit on 1Kevt of weighted data is closer to result of ML fit on 43Kevt of unweighted data 
00142   // than to 1Kevt of unweighted data, whereas the reference chi^2 fit with SumW2 error gives a result closer to
00143   // that of an unbinned ML fit to 1Kevt of unweighted data. 
00144 
00145   regResult(r_ml_unw10,"rf403_ml_unw10") ;
00146   regResult(r_ml_unw43,"rf403_ml_unw43") ;
00147   regResult(r_ml_wgt  ,"rf403_ml_wgt") ;
00148   regResult(r_chi2_wgt ,"rf403_ml_chi2") ;
00149   regPlot(frame,"rf403_plot1") ;
00150   
00151   delete binnedData ;
00152   delete data ;
00153   delete data2 ;
00154   delete data3 ;
00155 
00156   return kTRUE ;
00157   }
00158 } ;

Generated on Tue Jul 5 15:07:50 2011 for ROOT_528-00b_version by  doxygen 1.5.1