rf607_fitresult.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
00004 // 
00005 // Demonstration of options of the RooFitResult class
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 "RooChebychev.h"
00022 #include "RooFitResult.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 #include "TFile.h"
00027 #include "TStyle.h"
00028 #include "TH2.h"
00029 #include "TMatrixDSym.h"
00030 
00031 using namespace RooFit ;
00032 
00033 
00034 void rf607_fitresult()
00035 {
00036   // C r e a t e   p d f ,   d a t a
00037   // --------------------------------
00038 
00039   // Declare observable x
00040   RooRealVar x("x","x",0,10) ;
00041 
00042   // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their paramaters
00043   RooRealVar mean("mean","mean of gaussians",5,-10,10) ;
00044   RooRealVar sigma1("sigma1","width of gaussians",0.5,0.1,10) ;
00045   RooRealVar sigma2("sigma2","width of gaussians",1,0.1,10) ;
00046 
00047   RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;  
00048   RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;  
00049   
00050   // Build Chebychev polynomial p.d.f.  
00051   RooRealVar a0("a0","a0",0.5,0.,1.) ;
00052   RooRealVar a1("a1","a1",-0.2) ;
00053   RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
00054 
00055   // Sum the signal components into a composite signal p.d.f.
00056   RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
00057   RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
00058 
00059   // Sum the composite signal and background 
00060   RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
00061   RooAddPdf  model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
00062 
00063   // Generate 1000 events
00064   RooDataSet* data = model.generate(x,1000) ;
00065 
00066 
00067 
00068   // F i t   p d f   t o   d a t a ,   s a v e   f i t r e s u l t 
00069   // -------------------------------------------------------------
00070 
00071   // Perform fit and save result
00072   RooFitResult* r = model.fitTo(*data,Save()) ;
00073 
00074 
00075 
00076   // P r i n t   f i t   r e s u l t s 
00077   // ---------------------------------
00078 
00079   // Summary printing: Basic info plus final values of floating fit parameters
00080   r->Print() ;
00081 
00082   // Verbose printing: Basic info, values of constant paramaters, initial and
00083   // final values of floating parameters, global correlations
00084   r->Print("v") ;
00085 
00086 
00087 
00088   // V i s u a l i z e   c o r r e l a t i o n   m a t r i x
00089   // -------------------------------------------------------
00090 
00091   // Construct 2D color plot of correlation matrix
00092   gStyle->SetOptStat(0) ;
00093   gStyle->SetPalette(1) ;
00094   TH2* hcorr = r->correlationHist() ;
00095 
00096 
00097   // Visualize ellipse corresponding to single correlation matrix element
00098   RooPlot* frame = new RooPlot(sigma1,sig1frac,0.45,0.60,0.65,0.90) ;
00099   frame->SetTitle("Covariance between sigma1 and sig1frac") ;
00100   r->plotOn(frame,sigma1,sig1frac,"ME12ABHV") ;
00101 
00102 
00103 
00104   // A c c e s s   f i t   r e s u l t   i n f o r m a t i o n
00105   // ---------------------------------------------------------
00106 
00107   // Access basic information
00108   cout << "EDM = " << r->edm() << endl ;
00109   cout << "-log(L) at minimum = " << r->minNll() << endl ;
00110 
00111   // Access list of final fit parameter values
00112   cout << "final value of floating parameters" << endl ;
00113   r->floatParsFinal().Print("s") ;
00114 
00115   // Access correlation matrix elements
00116   cout << "correlation between sig1frac and a0 is  " << r->correlation(sig1frac,a0) << endl ;
00117   cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac","mean") << endl ;
00118 
00119   // Extract covariance and correlation matrix as TMatrixDSym
00120   const TMatrixDSym& cor = r->correlationMatrix() ;
00121   const TMatrixDSym& cov = r->covarianceMatrix() ;
00122 
00123   // Print correlation, covariance matrix
00124   cout << "correlation matrix" << endl ;
00125   cor.Print() ;
00126   cout << "covariance matrix" << endl ;
00127   cov.Print() ;
00128   
00129 
00130   // P e r s i s t   f i t   r e s u l t   i n   r o o t   f i l e 
00131   // -------------------------------------------------------------
00132   
00133   // Open new ROOT file save save result 
00134   TFile f("rf607_fitresult.root","RECREATE") ;
00135   r->Write("rf607") ;
00136   f.Close() ;
00137 
00138   // In a clean ROOT session retrieve the persisted fit result as follows:
00139   // RooFitResult* r = gDirectory->Get("rf607") ;
00140  
00141 
00142   TCanvas* c = new TCanvas("rf607_fitresult","rf607_fitresult",800,400) ;
00143   c->Divide(2) ;
00144   c->cd(1) ; gPad->SetLeftMargin(0.15) ; hcorr->GetYaxis()->SetTitleOffset(1.4) ; hcorr->Draw("colz") ;
00145   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00146 
00147 }

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