Logo ROOT  
Reference Guide
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 ///
10 /// \date July 2008
11 /// \author Wouter Verkerke
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "RooConstVar.h"
17 #include "RooAddPdf.h"
18 #include "RooChebychev.h"
19 #include "RooFitResult.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 #include "TFile.h"
24 #include "TStyle.h"
25 #include "TH2.h"
26 #include "TMatrixDSym.h"
27 
28 using namespace RooFit;
29 
30 void rf607_fitresult()
31 {
32  // C r e a t e p d f , d a t a
33  // --------------------------------
34 
35  // Declare observable x
36  RooRealVar x("x", "x", 0, 10);
37 
38  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
39  RooRealVar mean("mean", "mean of gaussians", 5, -10, 10);
40  RooRealVar sigma1("sigma1", "width of gaussians", 0.5, 0.1, 10);
41  RooRealVar sigma2("sigma2", "width of gaussians", 1, 0.1, 10);
42 
43  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
44  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
45 
46  // Build Chebychev polynomial pdf
47  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
48  RooRealVar a1("a1", "a1", -0.2);
49  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
50 
51  // Sum the signal components into a composite signal pdf
52  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
53  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
54 
55  // Sum the composite signal and background
56  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
57  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
58 
59  // Generate 1000 events
60  RooDataSet *data = model.generate(x, 1000);
61 
62  // 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
63  // -------------------------------------------------------------
64 
65  // Perform fit and save result
66  RooFitResult *r = model.fitTo(*data, Save());
67 
68  // P r i n t f i t r e s u l t s
69  // ---------------------------------
70 
71  // Summary printing: Basic info plus final values of floating fit parameters
72  r->Print();
73 
74  // Verbose printing: Basic info, values of constant parameters, initial and
75  // final values of floating parameters, global correlations
76  r->Print("v");
77 
78  // 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
79  // -------------------------------------------------------
80 
81  // Construct 2D color plot of correlation matrix
82  gStyle->SetOptStat(0);
83  TH2 *hcorr = r->correlationHist();
84 
85  // Visualize ellipse corresponding to single correlation matrix element
86  RooPlot *frame = new RooPlot(sigma1, sig1frac, 0.45, 0.60, 0.65, 0.90);
87  frame->SetTitle("Covariance between sigma1 and sig1frac");
88  r->plotOn(frame, sigma1, sig1frac, "ME12ABHV");
89 
90  // 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
91  // ---------------------------------------------------------
92 
93  // Access basic information
94  cout << "EDM = " << r->edm() << endl;
95  cout << "-log(L) at minimum = " << r->minNll() << endl;
96 
97  // Access list of final fit parameter values
98  cout << "final value of floating parameters" << endl;
99  r->floatParsFinal().Print("s");
100 
101  // Access correlation matrix elements
102  cout << "correlation between sig1frac and a0 is " << r->correlation(sig1frac, a0) << endl;
103  cout << "correlation between bkgfrac and mean is " << r->correlation("bkgfrac", "mean") << endl;
104 
105  // Extract covariance and correlation matrix as TMatrixDSym
106  const TMatrixDSym &cor = r->correlationMatrix();
107  const TMatrixDSym &cov = r->covarianceMatrix();
108 
109  // Print correlation, covariance matrix
110  cout << "correlation matrix" << endl;
111  cor.Print();
112  cout << "covariance matrix" << endl;
113  cov.Print();
114 
115  // 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
116  // -------------------------------------------------------------
117 
118  // Open new ROOT file save save result
119  TFile f("rf607_fitresult.root", "RECREATE");
120  r->Write("rf607");
121  f.Close();
122 
123  // In a clean ROOT session retrieve the persisted fit result as follows:
124  // RooFitResult* r = gDirectory->Get("rf607") ;
125 
126  TCanvas *c = new TCanvas("rf607_fitresult", "rf607_fitresult", 800, 400);
127  c->Divide(2);
128  c->cd(1);
129  gPad->SetLeftMargin(0.15);
130  hcorr->GetYaxis()->SetTitleOffset(1.4);
131  hcorr->Draw("colz");
132  c->cd(2);
133  gPad->SetLeftMargin(0.15);
134  frame->GetYaxis()->SetTitleOffset(1.6);
135  frame->Draw();
136 }
c
#define c(i)
Definition: RSha256.hxx:101
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
RooChebychev.h
RooAddPdf
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:32
f
#define f(i)
Definition: RSha256.hxx:104
TObject::Print
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:552
RooChebychev
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
r
ROOT::R::TRInterface & r
Definition: Object.C:4
RooArgList
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooGaussian.h
TStyle.h
x
Double_t x[n]
Definition: legend1.C:17
RooGaussian
Plain Gaussian p.d.f.
Definition: RooGaussian.h:24
RooAddPdf.h
TMatrixTSym< Double_t >
TCanvas.h
RooDataSet.h
TFile.h
RooFitResult
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
TMatrixDSym.h
rf607_fitresult
Definition: rf607_fitresult.py:1
gStyle
R__EXTERN TStyle * gStyle
Definition: TStyle.h:412
TH1::GetYaxis
TAxis * GetYaxis()
Definition: TH1.h:321
RooFit
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition: RooCFunction1Binding.h:29
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1258
TH2
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
RooPlot
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
RooRealVar.h
TH2.h
TStyle::SetOptStat
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:1589
TFile
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:54
RooFitResult.h
RooConstVar.h
RooPlot::SetTitle
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1237
TCanvas
The Canvas class.
Definition: TCanvas.h:23
TObject::Write
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:795
TAxis.h
gPad
#define gPad
Definition: TVirtualPad.h:287
RooDataSet
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
make_cnn_model.model
model
Definition: make_cnn_model.py:6
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:37
RooFit::Save
RooCmdArg Save(Bool_t flag=kTRUE)
Definition: RooGlobalFunc.cxx:190
RooArgSet
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:29
TH1::Draw
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:3050