Logo ROOT   6.07/09
Reference Guide
rf607_fitresult.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
5 ///
6 /// Demonstration of options of the RooFitResult class
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooAddPdf.h"
19 #include "RooChebychev.h"
20 #include "RooFitResult.h"
21 #include "TCanvas.h"
22 #include "TAxis.h"
23 #include "RooPlot.h"
24 #include "TFile.h"
25 #include "TStyle.h"
26 #include "TH2.h"
27 #include "TMatrixDSym.h"
28 
29 using namespace RooFit ;
30 
31 
32 void rf607_fitresult()
33 {
34  // C r e a t e p d f , d a t a
35  // --------------------------------
36 
37  // Declare observable x
38  RooRealVar x("x","x",0,10) ;
39 
40  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
41  RooRealVar mean("mean","mean of gaussians",5,-10,10) ;
42  RooRealVar sigma1("sigma1","width of gaussians",0.5,0.1,10) ;
43  RooRealVar sigma2("sigma2","width of gaussians",1,0.1,10) ;
44 
45  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
46  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
47 
48  // Build Chebychev polynomial p.d.f.
49  RooRealVar a0("a0","a0",0.5,0.,1.) ;
50  RooRealVar a1("a1","a1",-0.2) ;
51  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
52 
53  // Sum the signal components into a composite signal p.d.f.
54  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
55  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
56 
57  // Sum the composite signal and background
58  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
59  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
60 
61  // Generate 1000 events
62  RooDataSet* data = model.generate(x,1000) ;
63 
64 
65 
66  // 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
67  // -------------------------------------------------------------
68 
69  // Perform fit and save result
70  RooFitResult* r = model.fitTo(*data,Save()) ;
71 
72 
73 
74  // P r i n t f i t r e s u l t s
75  // ---------------------------------
76 
77  // Summary printing: Basic info plus final values of floating fit parameters
78  r->Print() ;
79 
80  // Verbose printing: Basic info, values of constant parameters, initial and
81  // final values of floating parameters, global correlations
82  r->Print("v") ;
83 
84 
85 
86  // 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
87  // -------------------------------------------------------
88 
89  // Construct 2D color plot of correlation matrix
90  gStyle->SetOptStat(0) ;
91  TH2* hcorr = r->correlationHist() ;
92 
93 
94  // Visualize ellipse corresponding to single correlation matrix element
95  RooPlot* frame = new RooPlot(sigma1,sig1frac,0.45,0.60,0.65,0.90) ;
96  frame->SetTitle("Covariance between sigma1 and sig1frac") ;
97  r->plotOn(frame,sigma1,sig1frac,"ME12ABHV") ;
98 
99 
100 
101  // 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
102  // ---------------------------------------------------------
103 
104  // Access basic information
105  cout << "EDM = " << r->edm() << endl ;
106  cout << "-log(L) at minimum = " << r->minNll() << endl ;
107 
108  // Access list of final fit parameter values
109  cout << "final value of floating parameters" << endl ;
110  r->floatParsFinal().Print("s") ;
111 
112  // Access correlation matrix elements
113  cout << "correlation between sig1frac and a0 is " << r->correlation(sig1frac,a0) << endl ;
114  cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac","mean") << endl ;
115 
116  // Extract covariance and correlation matrix as TMatrixDSym
117  const TMatrixDSym& cor = r->correlationMatrix() ;
118  const TMatrixDSym& cov = r->covarianceMatrix() ;
119 
120  // Print correlation, covariance matrix
121  cout << "correlation matrix" << endl ;
122  cor.Print() ;
123  cout << "covariance matrix" << endl ;
124  cov.Print() ;
125 
126 
127  // 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
128  // -------------------------------------------------------------
129 
130  // Open new ROOT file save save result
131  TFile f("rf607_fitresult.root","RECREATE") ;
132  r->Write("rf607") ;
133  f.Close() ;
134 
135  // In a clean ROOT session retrieve the persisted fit result as follows:
136  // RooFitResult* r = gDirectory->Get("rf607") ;
137 
138 
139  TCanvas* c = new TCanvas("rf607_fitresult","rf607_fitresult",800,400) ;
140  c->Divide(2) ;
141  c->cd(1) ; gPad->SetLeftMargin(0.15) ; hcorr->GetYaxis()->SetTitleOffset(1.4) ; hcorr->Draw("colz") ;
142  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
143 
144 }
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:262
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:830
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
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
return c
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:50
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
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 &#39;title&#39;.
Definition: RooPlot.cxx:1099
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
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.
TRandom2 r(17)
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:2853
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:325
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
The Canvas class.
Definition: TCanvas.h:41
double f(double x)
Double_t minNll() const
Definition: RooFitResult.h:97
RooCmdArg Save(Bool_t flag=kTRUE)
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1089
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:1257
#define gPad
Definition: TVirtualPad.h:289
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const TMatrixDSym & correlationMatrix() const
Return correlation matrix ;.
Double_t edm() const
Definition: RooFitResult.h:93
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559