ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf607_fitresult.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
4 ///
5 /// Demonstration of options of the RooFitResult class
6 ///
7 ///
8 ///
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooConstVar.h"
20 #include "RooAddPdf.h"
21 #include "RooChebychev.h"
22 #include "RooFitResult.h"
23 #include "TCanvas.h"
24 #include "TAxis.h"
25 #include "RooPlot.h"
26 #include "TFile.h"
27 #include "TStyle.h"
28 #include "TH2.h"
29 #include "TMatrixDSym.h"
30 
31 using namespace RooFit ;
32 
33 
34 void rf607_fitresult()
35 {
36  // C r e a t e p d f , d a t a
37  // --------------------------------
38 
39  // Declare observable x
40  RooRealVar x("x","x",0,10) ;
41 
42  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
43  RooRealVar mean("mean","mean of gaussians",5,-10,10) ;
44  RooRealVar sigma1("sigma1","width of gaussians",0.5,0.1,10) ;
45  RooRealVar sigma2("sigma2","width of gaussians",1,0.1,10) ;
46 
47  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
48  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
49 
50  // Build Chebychev polynomial p.d.f.
51  RooRealVar a0("a0","a0",0.5,0.,1.) ;
52  RooRealVar a1("a1","a1",-0.2) ;
53  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
54 
55  // Sum the signal components into a composite signal p.d.f.
56  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
57  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
58 
59  // Sum the composite signal and background
60  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
61  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
62 
63  // Generate 1000 events
64  RooDataSet* data = model.generate(x,1000) ;
65 
66 
67 
68  // 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
69  // -------------------------------------------------------------
70 
71  // Perform fit and save result
72  RooFitResult* r = model.fitTo(*data,Save()) ;
73 
74 
75 
76  // P r i n t f i t r e s u l t s
77  // ---------------------------------
78 
79  // Summary printing: Basic info plus final values of floating fit parameters
80  r->Print() ;
81 
82  // Verbose printing: Basic info, values of constant parameters, initial and
83  // final values of floating parameters, global correlations
84  r->Print("v") ;
85 
86 
87 
88  // 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
89  // -------------------------------------------------------
90 
91  // Construct 2D color plot of correlation matrix
92  gStyle->SetOptStat(0) ;
93  gStyle->SetPalette(1) ;
94  TH2* hcorr = r->correlationHist() ;
95 
96 
97  // Visualize ellipse corresponding to single correlation matrix element
98  RooPlot* frame = new RooPlot(sigma1,sig1frac,0.45,0.60,0.65,0.90) ;
99  frame->SetTitle("Covariance between sigma1 and sig1frac") ;
100  r->plotOn(frame,sigma1,sig1frac,"ME12ABHV") ;
101 
102 
103 
104  // 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
105  // ---------------------------------------------------------
106 
107  // Access basic information
108  cout << "EDM = " << r->edm() << endl ;
109  cout << "-log(L) at minimum = " << r->minNll() << endl ;
110 
111  // Access list of final fit parameter values
112  cout << "final value of floating parameters" << endl ;
113  r->floatParsFinal().Print("s") ;
114 
115  // Access correlation matrix elements
116  cout << "correlation between sig1frac and a0 is " << r->correlation(sig1frac,a0) << endl ;
117  cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac","mean") << endl ;
118 
119  // Extract covariance and correlation matrix as TMatrixDSym
120  const TMatrixDSym& cor = r->correlationMatrix() ;
121  const TMatrixDSym& cov = r->covarianceMatrix() ;
122 
123  // Print correlation, covariance matrix
124  cout << "correlation matrix" << endl ;
125  cor.Print() ;
126  cout << "covariance matrix" << endl ;
127  cov.Print() ;
128 
129 
130  // 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
131  // -------------------------------------------------------------
132 
133  // Open new ROOT file save save result
134  TFile f("rf607_fitresult.root","RECREATE") ;
135  r->Write("rf607") ;
136  f.Close() ;
137 
138  // In a clean ROOT session retrieve the persisted fit result as follows:
139  // RooFitResult* r = gDirectory->Get("rf607") ;
140 
141 
142  TCanvas* c = new TCanvas("rf607_fitresult","rf607_fitresult",800,400) ;
143  c->Divide(2) ;
144  c->cd(1) ; gPad->SetLeftMargin(0.15) ; hcorr->GetYaxis()->SetTitleOffset(1.4) ; hcorr->Draw("colz") ;
145  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
146 
147 }
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:65
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:245
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
RooPlot * plotOn(RooPlot *frame, const RooAbsArg &par1, const RooAbsArg &par2, const char *options="ME") const
Definition: RooFitResult.h:142
const RooArgList & floatParsFinal() const
Definition: RooFitResult.h:109
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
TH2 * correlationHist(const char *name="correlation_matrix") const
Return TH2D of correlation matrix.
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1099
TFile * f
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
friend class RooArgSet
Definition: RooAbsArg.h:469
Double_t correlation(const RooAbsArg &par1, const RooAbsArg &par2) const
Definition: RooFitResult.h:116
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
void Print(Option_t *name="") const
Print the matrix as a table of elements.
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
TAxis * GetYaxis()
Definition: TH1.h:320
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
Double_t minNll() const
Definition: RooFitResult.h:97
RooCmdArg Save(Bool_t flag=kTRUE)
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1252
#define gPad
Definition: TVirtualPad.h:288
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const TMatrixDSym & correlationMatrix() const
Return correlation matrix ;.
Double_t edm() const
Definition: RooFitResult.h:93
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition: TStyle.cxx:1445
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:297
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559