rf605_profilell.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605
00004 // 
00005 // Working with the profile likelihood estimator
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 "RooGaussian.h"
00019 #include "RooAddPdf.h"
00020 #include "RooNLLVar.h"
00021 #include "RooProfileLL.h"
00022 #include "RooMinuit.h"
00023 #include "TCanvas.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026 
00027 
00028 class TestBasic605 : public RooFitTestUnit
00029 {
00030 public: 
00031   TestBasic605(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Profile Likelihood operator",refFile,writeRef,verbose) {} ;
00032   Bool_t testCode() {
00033 
00034   // C r e a t e   m o d e l   a n d   d a t a s e t 
00035   // -----------------------------------------------
00036 
00037   // Observable
00038   RooRealVar x("x","x",-20,20) ;
00039 
00040   // Model (intentional strong correlations)
00041   RooRealVar mean("mean","mean of g1 and g2",0,-10,10) ;
00042   RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
00043   RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
00044 
00045   RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
00046   RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
00047 
00048   RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
00049   RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
00050 
00051   // Generate 1000 events
00052   RooDataSet* data = model.generate(x,1000) ;
00053   
00054 
00055 
00056   // C o n s t r u c t   p l a i n   l i k e l i h o o d
00057   // ---------------------------------------------------
00058 
00059   // Construct unbinned likelihood
00060   RooNLLVar nll("nll","nll",model,*data) ;
00061 
00062   // Minimize likelihood w.r.t all parameters before making plots
00063   RooMinuit(nll).migrad() ;
00064 
00065   // Plot likelihood scan frac 
00066   RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
00067   nll.plotOn(frame1,ShiftToZero()) ;
00068 
00069   // Plot likelihood scan in sigma_g2
00070   RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
00071   nll.plotOn(frame2,ShiftToZero()) ;
00072 
00073 
00074 
00075   // C o n s t r u c t   p r o f i l e   l i k e l i h o o d   i n   f r a c
00076   // -----------------------------------------------------------------------
00077 
00078   // The profile likelihood estimator on nll for frac will minimize nll w.r.t
00079   // all floating parameters except frac for each evaluation
00080   RooProfileLL pll_frac("pll_frac","pll_frac",nll,frac) ;
00081 
00082   // Plot the profile likelihood in frac
00083   pll_frac.plotOn(frame1,LineColor(kRed)) ;
00084 
00085   // Adjust frame maximum for visual clarity
00086   frame1->SetMinimum(0) ;
00087   frame1->SetMaximum(3) ;
00088 
00089 
00090 
00091   // C o n s t r u c t   p r o f i l e   l i k e l i h o o d   i n   s i g m a _ g 2 
00092   // -------------------------------------------------------------------------------
00093 
00094   // The profile likelihood estimator on nll for sigma_g2 will minimize nll
00095   // w.r.t all floating parameters except sigma_g2 for each evaluation
00096   RooProfileLL pll_sigmag2("pll_sigmag2","pll_sigmag2",nll,sigma_g2) ;
00097 
00098   // Plot the profile likelihood in sigma_g2
00099   pll_sigmag2.plotOn(frame2,LineColor(kRed)) ;
00100 
00101   // Adjust frame maximum for visual clarity
00102   frame2->SetMinimum(0) ;
00103   frame2->SetMaximum(3) ;
00104 
00105 
00106   regPlot(frame1,"rf605_plot1") ;
00107   regPlot(frame2,"rf605_plot2") ;
00108 
00109   delete data ;
00110 
00111   return kTRUE ;
00112   }
00113 } ;
00114 

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