Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf608_fitresultaspdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// Likelihood and minimization: representing the parabolic approximation of the fit as a multi-variate Gaussian on the
5 /// parameters of the fitted p.d.f.
6 ///
7 /// \macro_image
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 #include "RooRealVar.h"
13 #include "RooDataSet.h"
14 #include "RooGaussian.h"
15 #include "RooConstVar.h"
16 #include "RooAddPdf.h"
17 #include "RooChebychev.h"
18 #include "RooFitResult.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "RooPlot.h"
22 #include "TFile.h"
23 #include "TStyle.h"
24 #include "TH2.h"
25 #include "TH3.h"
26 
27 using namespace RooFit;
28 
29 void rf608_fitresultaspdf()
30 {
31  // C r e a t e m o d e l a n d d a t a s e t
32  // -----------------------------------------------
33 
34  // Observable
35  RooRealVar x("x", "x", -20, 20);
36 
37  // Model (intentional strong correlations)
38  RooRealVar mean("mean", "mean of g1 and g2", 0, -1, 1);
39  RooRealVar sigma_g1("sigma_g1", "width of g1", 2);
40  RooGaussian g1("g1", "g1", x, mean, sigma_g1);
41 
42  RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 5.0);
43  RooGaussian g2("g2", "g2", x, mean, sigma_g2);
44 
45  RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
46  RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
47 
48  // Generate 1000 events
49  RooDataSet *data = model.generate(x, 1000);
50 
51  // F i t m o d e l t o d a t a
52  // ----------------------------------
53 
54  RooFitResult *r = model.fitTo(*data, Save());
55 
56  // C r e a t e M V G a u s s i a n p d f o f f i t t e d p a r a m e t e r s
57  // ------------------------------------------------------------------------------------
58 
59  RooAbsPdf *parabPdf = r->createHessePdf(RooArgSet(frac, mean, sigma_g2));
60 
61  // S o m e e x e c e r c i s e s w i t h t h e p a r a m e t e r p d f
62  // -----------------------------------------------------------------------------
63 
64  // Generate 100K points in the parameter space, sampled from the MVGaussian p.d.f.
65  RooDataSet *d = parabPdf->generate(RooArgSet(mean, sigma_g2, frac), 100000);
66 
67  // Sample a 3-D histogram of the p.d.f. to be visualized as an error ellipsoid using the GLISO draw option
68  TH3 *hh_3d = (TH3 *)parabPdf->createHistogram("mean,sigma_g2,frac", 25, 25, 25);
69  hh_3d->SetFillColor(kBlue);
70 
71  // Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s
72  // The integrations corresponding to these projections are performed analytically
73  // by the MV Gaussian p.d.f.
74  RooAbsPdf *pdf_sigmag2_frac = parabPdf->createProjection(mean);
75  RooAbsPdf *pdf_mean_frac = parabPdf->createProjection(sigma_g2);
76  RooAbsPdf *pdf_mean_sigmag2 = parabPdf->createProjection(frac);
77 
78  // Make 2D plots of the 3 two-dimensional p.d.f. projections
79  TH2 *hh_sigmag2_frac = (TH2 *)pdf_sigmag2_frac->createHistogram("sigma_g2,frac", 50, 50);
80  TH2 *hh_mean_frac = (TH2 *)pdf_mean_frac->createHistogram("mean,frac", 50, 50);
81  TH2 *hh_mean_sigmag2 = (TH2 *)pdf_mean_sigmag2->createHistogram("mean,sigma_g2", 50, 50);
82  hh_mean_frac->SetLineColor(kBlue);
83  hh_sigmag2_frac->SetLineColor(kBlue);
84  hh_mean_sigmag2->SetLineColor(kBlue);
85 
86  // Draw the 'sigar'
87  new TCanvas("rf608_fitresultaspdf_1", "rf608_fitresultaspdf_1", 600, 600);
88  hh_3d->Draw("iso");
89 
90  // Draw the 2D projections of the 3D p.d.f.
91  TCanvas *c2 = new TCanvas("rf608_fitresultaspdf_2", "rf608_fitresultaspdf_2", 900, 600);
92  c2->Divide(3, 2);
93  c2->cd(1);
94  gPad->SetLeftMargin(0.15);
95  hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4);
96  hh_mean_sigmag2->Draw("surf3");
97  c2->cd(2);
98  gPad->SetLeftMargin(0.15);
99  hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4);
100  hh_sigmag2_frac->Draw("surf3");
101  c2->cd(3);
102  gPad->SetLeftMargin(0.15);
103  hh_mean_frac->GetZaxis()->SetTitleOffset(1.4);
104  hh_mean_frac->Draw("surf3");
105 
106  // Draw the distributions of parameter points sampled from the p.d.f.
107  TH1 *tmp1 = d->createHistogram("mean,sigma_g2", 50, 50);
108  TH1 *tmp2 = d->createHistogram("sigma_g2,frac", 50, 50);
109  TH1 *tmp3 = d->createHistogram("mean,frac", 50, 50);
110 
111  c2->cd(4);
112  gPad->SetLeftMargin(0.15);
113  tmp1->GetZaxis()->SetTitleOffset(1.4);
114  tmp1->Draw("lego3");
115  c2->cd(5);
116  gPad->SetLeftMargin(0.15);
117  tmp2->GetZaxis()->SetTitleOffset(1.4);
118  tmp2->Draw("lego3");
119  c2->cd(6);
120  gPad->SetLeftMargin(0.15);
121  tmp3->GetZaxis()->SetTitleOffset(1.4);
122  tmp3->Draw("lego3");
123 }