rf605_profilell.C

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 "RooConstVar.h"
00020 #include "RooAddPdf.h"
00021 #include "RooMinuit.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026 
00027 
00028 void rf605_profilell()
00029 {
00030   // C r e a t e   m o d e l   a n d   d a t a s e t 
00031   // -----------------------------------------------
00032 
00033   // Observable
00034   RooRealVar x("x","x",-20,20) ;
00035 
00036   // Model (intentional strong correlations)
00037   RooRealVar mean("mean","mean of g1 and g2",0,-10,10) ;
00038   RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
00039   RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
00040 
00041   RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
00042   RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
00043 
00044   RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
00045   RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
00046 
00047   // Generate 1000 events
00048   RooDataSet* data = model.generate(x,1000) ;
00049   
00050 
00051 
00052   // C o n s t r u c t   p l a i n   l i k e l i h o o d
00053   // ---------------------------------------------------
00054 
00055   // Construct unbinned likelihood
00056   RooAbsReal* nll = model.createNLL(*data,NumCPU(2)) ;
00057 
00058   // Minimize likelihood w.r.t all parameters before making plots
00059   RooMinuit(*nll).migrad() ;
00060 
00061   // Plot likelihood scan frac 
00062   RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
00063   nll->plotOn(frame1,ShiftToZero()) ;
00064 
00065   // Plot likelihood scan in sigma_g2
00066   RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
00067   nll->plotOn(frame2,ShiftToZero()) ;
00068 
00069 
00070 
00071   // 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
00072   // -----------------------------------------------------------------------
00073 
00074   // The profile likelihood estimator on nll for frac will minimize nll w.r.t
00075   // all floating parameters except frac for each evaluation
00076 
00077   RooAbsReal* pll_frac = nll->createProfile(frac) ;
00078 
00079   // Plot the profile likelihood in frac
00080   pll_frac->plotOn(frame1,LineColor(kRed)) ;
00081 
00082   // Adjust frame maximum for visual clarity
00083   frame1->SetMinimum(0) ;
00084   frame1->SetMaximum(3) ;
00085 
00086 
00087 
00088   // 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 
00089   // -------------------------------------------------------------------------------
00090 
00091   // The profile likelihood estimator on nll for sigma_g2 will minimize nll
00092   // w.r.t all floating parameters except sigma_g2 for each evaluation
00093   RooAbsReal* pll_sigmag2 = nll->createProfile(sigma_g2) ;
00094 
00095   // Plot the profile likelihood in sigma_g2
00096   pll_sigmag2->plotOn(frame2,LineColor(kRed)) ;
00097 
00098   // Adjust frame maximum for visual clarity
00099   frame2->SetMinimum(0) ;
00100   frame2->SetMaximum(3) ;
00101 
00102 
00103 
00104   // Make canvas and draw RooPlots
00105   TCanvas *c = new TCanvas("rf605_profilell","rf605_profilell",800, 400);
00106   c->Divide(2);
00107   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
00108   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
00109 
00110   delete pll_frac ;
00111   delete pll_sigmag2 ;
00112   delete nll ;
00113 }
00114 

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