Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf607_fitresult.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// Likelihood and minimization: demonstration of options of the RooFitResult class
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooAddPdf.h"
16 #include "RooChebychev.h"
17 #include "RooFitResult.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 #include "TFile.h"
22 #include "TStyle.h"
23 #include "TH2.h"
24 #include "TMatrixDSym.h"
25 
26 using namespace RooFit;
27 
28 void rf607_fitresult()
29 {
30  // C r e a t e p d f , d a t a
31  // --------------------------------
32 
33  // Declare observable x
34  RooRealVar x("x", "x", 0, 10);
35 
36  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
37  RooRealVar mean("mean", "mean of gaussians", 5, -10, 10);
38  RooRealVar sigma1("sigma1", "width of gaussians", 0.5, 0.1, 10);
39  RooRealVar sigma2("sigma2", "width of gaussians", 1, 0.1, 10);
40 
41  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
42  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
43 
44  // Build Chebychev polynomial p.d.f.
45  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
46  RooRealVar a1("a1", "a1", -0.2);
47  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
48 
49  // Sum the signal components into a composite signal p.d.f.
50  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
51  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
52 
53  // Sum the composite signal and background
54  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
55  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
56 
57  // Generate 1000 events
58  RooDataSet *data = model.generate(x, 1000);
59 
60  // 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
61  // -------------------------------------------------------------
62 
63  // Perform fit and save result
64  RooFitResult *r = model.fitTo(*data, Save());
65 
66  // P r i n t f i t r e s u l t s
67  // ---------------------------------
68 
69  // Summary printing: Basic info plus final values of floating fit parameters
70  r->Print();
71 
72  // Verbose printing: Basic info, values of constant parameters, initial and
73  // final values of floating parameters, global correlations
74  r->Print("v");
75 
76  // 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
77  // -------------------------------------------------------
78 
79  // Construct 2D color plot of correlation matrix
80  gStyle->SetOptStat(0);
81  TH2 *hcorr = r->correlationHist();
82 
83  // Visualize ellipse corresponding to single correlation matrix element
84  RooPlot *frame = new RooPlot(sigma1, sig1frac, 0.45, 0.60, 0.65, 0.90);
85  frame->SetTitle("Covariance between sigma1 and sig1frac");
86  r->plotOn(frame, sigma1, sig1frac, "ME12ABHV");
87 
88  // 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
89  // ---------------------------------------------------------
90 
91  // Access basic information
92  cout << "EDM = " << r->edm() << endl;
93  cout << "-log(L) at minimum = " << r->minNll() << endl;
94 
95  // Access list of final fit parameter values
96  cout << "final value of floating parameters" << endl;
97  r->floatParsFinal().Print("s");
98 
99  // Access correlation matrix elements
100  cout << "correlation between sig1frac and a0 is " << r->correlation(sig1frac, a0) << endl;
101  cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac", "mean") << endl;
102 
103  // Extract covariance and correlation matrix as TMatrixDSym
104  const TMatrixDSym &cor = r->correlationMatrix();
105  const TMatrixDSym &cov = r->covarianceMatrix();
106 
107  // Print correlation, covariance matrix
108  cout << "correlation matrix" << endl;
109  cor.Print();
110  cout << "covariance matrix" << endl;
111  cov.Print();
112 
113  // 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
114  // -------------------------------------------------------------
115 
116  // Open new ROOT file save save result
117  TFile f("rf607_fitresult.root", "RECREATE");
118  r->Write("rf607");
119  f.Close();
120 
121  // In a clean ROOT session retrieve the persisted fit result as follows:
122  // RooFitResult* r = gDirectory->Get("rf607") ;
123 
124  TCanvas *c = new TCanvas("rf607_fitresult", "rf607_fitresult", 800, 400);
125  c->Divide(2);
126  c->cd(1);
127  gPad->SetLeftMargin(0.15);
128  hcorr->GetYaxis()->SetTitleOffset(1.4);
129  hcorr->Draw("colz");
130  c->cd(2);
131  gPad->SetLeftMargin(0.15);
132  frame->GetYaxis()->SetTitleOffset(1.6);
133  frame->Draw();
134 }