rf606_nllerrorhandling.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #606
00004 // 
00005 // Understanding and customizing error handling in likelihood evaluations
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 "RooArgusBG.h"
00019 #include "RooNLLVar.h"
00020 #include "TCanvas.h"
00021 #include "RooPlot.h"
00022 using namespace RooFit ;
00023 
00024 
00025 class TestBasic606 : public RooFitTestUnit
00026 {
00027 public: 
00028   TestBasic606(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("NLL error handling",refFile,writeRef,verbose) {} ;
00029   Bool_t testCode() {
00030 
00031   // C r e a t e   m o d e l  a n d   d a t a s e t 
00032   // ----------------------------------------------
00033 
00034   // Observable
00035   RooRealVar m("m","m",5.20,5.30) ;
00036 
00037   // Parameters
00038   RooRealVar m0("m0","m0",5.291,5.20,5.30) ;
00039   RooRealVar k("k","k",-30,-50,-10) ;
00040 
00041   // Pdf
00042   RooArgusBG argus("argus","argus",m,m0,k) ;
00043 
00044   // Sample 1000 events in m from argus
00045   RooDataSet* data = argus.generate(m,1000) ;
00046 
00047 
00048 
00049   // P l o t   m o d e l   a n d   d a t a 
00050   // --------------------------------------
00051 
00052   RooPlot* frame1 = m.frame(Bins(40),Title("Argus model and data")) ;
00053   data->plotOn(frame1) ;
00054   argus.plotOn(frame1) ;
00055 
00056 
00057 
00058   // F i t   m o d e l   t o   d a t a 
00059   // ---------------------------------
00060 
00061   argus.fitTo(*data,PrintEvalErrors(10),Warnings(kFALSE)) ;    
00062   m0.setError(0.1) ;
00063   argus.fitTo(*data,PrintEvalErrors(0),EvalErrorWall(kFALSE),Warnings(kFALSE)) ;    
00064 
00065 
00066 
00067   // P l o t   l i k e l i h o o d   a s   f u n c t i o n   o f   m 0 
00068   // ------------------------------------------------------------------
00069 
00070   // Construct likelihood function of model and data
00071   RooNLLVar nll("nll","nll",argus,*data) ;
00072 
00073   // Plot likelihood in m0 in range that includes problematic values
00074   // In this configuration the number of errors per likelihood point 
00075   // evaluated for the curve is shown. A positive number in PrintEvalErrors(N)
00076   // will show details for up to N events. By default the values for likelihood
00077   // evaluations with errors are shown normally (unlike fitting), but the shape
00078   // of the curve can be erratic in these regions.
00079 
00080   RooPlot* frame2 = m0.frame(Range(5.288,5.293),Title("-log(L) scan vs m0")) ;
00081   nll.plotOn(frame2,PrintEvalErrors(0),ShiftToZero(),LineColor(kRed),Precision(1e-4)) ; 
00082 
00083 
00084   // Plot likelihood in m0 in range that includes problematic values
00085   // In this configuration no messages are printed for likelihood evaluation errors,
00086   // but if an likelihood value evaluates with error, the corresponding value
00087   // on the curve will be set to the value given in EvalErrorValue().
00088 
00089   RooPlot* frame3 = m0.frame(Range(5.288,5.293),Title("-log(L) scan vs m0, problematic regions masked")) ;
00090   nll.plotOn(frame3,PrintEvalErrors(-1),ShiftToZero(),EvalErrorValue(nll.getVal()+10),LineColor(kRed)) ; 
00091 
00092 
00093   regPlot(frame1,"rf606_plot1") ;
00094   regPlot(frame2,"rf606_plot2") ;
00095   regPlot(frame3,"rf606_plot3") ;
00096 
00097   delete data ;
00098   return kTRUE ;
00099 
00100   }
00101 } ;

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