rf601_intminuit.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #601
00004 // 
00005 // Interactive minimization with MINUIT
00006 //
00007 //
00008 // 07/2008 - Wouter Verkerke 
00009 //
00010 /////////////////////////////////////////////////////////////////////////
00011 
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooGaussian.h"
00018 #include "RooConstVar.h"
00019 #include "RooProdPdf.h"
00020 #include "RooAddPdf.h"
00021 #include "RooMinuit.h"
00022 #include "RooFitResult.h"
00023 #include "RooPlot.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "TH1.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf601_intminuit()
00031 {
00032   // S e t u p   p d f   a n d   l i k e l i h o o d 
00033   // -----------------------------------------------
00034 
00035   // Observable
00036   RooRealVar x("x","x",-20,20) ;
00037 
00038   // Model (intentional strong correlations)
00039   RooRealVar mean("mean","mean of g1 and g2",0) ;
00040   RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
00041   RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
00042 
00043   RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
00044   RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
00045 
00046   RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
00047   RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
00048 
00049   // Generate 1000 events
00050   RooDataSet* data = model.generate(x,1000) ;
00051   
00052   // Construct unbinned likelihood of model w.r.t. data
00053   RooAbsReal* nll = model.createNLL(*data) ;
00054 
00055   // I n t e r a c t i v e   m i n i m i z a t i o n ,   e r r o r   a n a l y s i s
00056   // -------------------------------------------------------------------------------
00057 
00058   // Create MINUIT interface object
00059   RooMinuit m(*nll) ;
00060 
00061   // Activate verbose logging of MINUIT parameter space stepping
00062   m.setVerbose(kTRUE) ;
00063 
00064   // Call MIGRAD to minimize the likelihood
00065   m.migrad() ;
00066 
00067   // Print values of all paramaters, that reflect values (and error estimates)
00068   // that are back propagated from MINUIT
00069   model.getParameters(x)->Print("s") ;
00070 
00071   // Disable verbose logging
00072   m.setVerbose(kFALSE) ;
00073 
00074   // Run HESSE to calculate errors from d2L/dp2
00075   m.hesse() ;
00076 
00077   // Print value (and error) of sigma_g2 parameter, that reflects
00078   // value and error back propagated from MINUIT
00079   sigma_g2.Print() ;
00080 
00081   // Run MINOS on sigma_g2 parameter only
00082   m.minos(sigma_g2) ;
00083 
00084   // Print value (and error) of sigma_g2 parameter, that reflects
00085   // value and error back propagated from MINUIT
00086   sigma_g2.Print() ;
00087 
00088 
00089 
00090   // S a v i n g   r e s u l t s ,   c o n t o u r   p l o t s 
00091   // ---------------------------------------------------------
00092 
00093   // Save a snapshot of the fit result. This object contains the initial
00094   // fit parameters, the final fit parameters, the complete correlation
00095   // matrix, the EDM, the minimized FCN , the last MINUIT status code and
00096   // the number of times the RooFit function object has indicated evaluation
00097   // problems (e.g. zero probabilities during likelihood evaluation)
00098   RooFitResult* r = m.save() ;
00099 
00100   // Make contour plot of mx vs sx at 1,2,3 sigma
00101   RooPlot* frame = m.contour(frac,sigma_g2,1,2,3) ;
00102   frame->SetTitle("RooMinuit contour plot") ;
00103 
00104   // Print the fit result snapshot
00105   r->Print("v") ;
00106 
00107 
00108 
00109   // C h a n g e   p a r a m e t e r   v a l u e s ,   f l o a t i n g
00110   // -----------------------------------------------------------------
00111 
00112   // At any moment you can manually change the value of a (constant)
00113   // parameter
00114   mean = 0.3 ;
00115   
00116   // Rerun MIGRAD,HESSE
00117   m.migrad() ;
00118   m.hesse() ;
00119   frac.Print() ;
00120 
00121   // Now fix sigma_g2
00122   sigma_g2.setConstant(kTRUE) ;
00123 
00124   // Rerun MIGRAD,HESSE
00125   m.migrad() ;
00126   m.hesse() ;
00127   frac.Print() ;
00128 
00129 
00130 
00131   new TCanvas("rf601_intminuit","rf601_intminuit",600,600) ;
00132   gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00133 
00134 }

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