Logo ROOT  
Reference Guide
rf608_fitresultaspdf.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4///
5/// Likelihood and minimization: representing the parabolic approximation of the fit as a multi-variate Gaussian on the
6/// parameters of the fitted p.d.f.
7///
8/// \macro_image
9/// \macro_output
10/// \macro_code
11///
12/// \date 07/2008
13/// \author Wouter Verkerke
14
15#include "RooRealVar.h"
16#include "RooDataSet.h"
17#include "RooGaussian.h"
18#include "RooConstVar.h"
19#include "RooAddPdf.h"
20#include "RooChebychev.h"
21#include "RooFitResult.h"
22#include "TCanvas.h"
23#include "TAxis.h"
24#include "RooPlot.h"
25#include "TFile.h"
26#include "TStyle.h"
27#include "TH2.h"
28#include "TH3.h"
29
30using namespace RooFit;
31
33{
34 // C r e a t e m o d e l a n d d a t a s e t
35 // -----------------------------------------------
36
37 // Observable
38 RooRealVar x("x", "x", -20, 20);
39
40 // Model (intentional strong correlations)
41 RooRealVar mean("mean", "mean of g1 and g2", 0, -1, 1);
42 RooRealVar sigma_g1("sigma_g1", "width of g1", 2);
43 RooGaussian g1("g1", "g1", x, mean, sigma_g1);
44
45 RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 5.0);
46 RooGaussian g2("g2", "g2", x, mean, sigma_g2);
47
48 RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
49 RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
50
51 // Generate 1000 events
52 RooDataSet *data = model.generate(x, 1000);
53
54 // F i t m o d e l t o d a t a
55 // ----------------------------------
56
57 RooFitResult *r = model.fitTo(*data, Save());
58
59 // 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
60 // ------------------------------------------------------------------------------------
61
62 RooAbsPdf *parabPdf = r->createHessePdf(RooArgSet(frac, mean, sigma_g2));
63
64 // 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
65 // -----------------------------------------------------------------------------
66
67 // Generate 100K points in the parameter space, sampled from the MVGaussian p.d.f.
68 RooDataSet *d = parabPdf->generate(RooArgSet(mean, sigma_g2, frac), 100000);
69
70 // Sample a 3-D histogram of the p.d.f. to be visualized as an error ellipsoid using the GLISO draw option
71 TH3 *hh_3d = (TH3 *)parabPdf->createHistogram("mean,sigma_g2,frac", 25, 25, 25);
72 hh_3d->SetFillColor(kBlue);
73
74 // Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s
75 // The integrations corresponding to these projections are performed analytically
76 // by the MV Gaussian p.d.f.
77 RooAbsPdf *pdf_sigmag2_frac = parabPdf->createProjection(mean);
78 RooAbsPdf *pdf_mean_frac = parabPdf->createProjection(sigma_g2);
79 RooAbsPdf *pdf_mean_sigmag2 = parabPdf->createProjection(frac);
80
81 // Make 2D plots of the 3 two-dimensional p.d.f. projections
82 TH2 *hh_sigmag2_frac = (TH2 *)pdf_sigmag2_frac->createHistogram("sigma_g2,frac", 50, 50);
83 TH2 *hh_mean_frac = (TH2 *)pdf_mean_frac->createHistogram("mean,frac", 50, 50);
84 TH2 *hh_mean_sigmag2 = (TH2 *)pdf_mean_sigmag2->createHistogram("mean,sigma_g2", 50, 50);
85 hh_mean_frac->SetLineColor(kBlue);
86 hh_sigmag2_frac->SetLineColor(kBlue);
87 hh_mean_sigmag2->SetLineColor(kBlue);
88
89 // Draw the 'sigar'
90 new TCanvas("rf608_fitresultaspdf_1", "rf608_fitresultaspdf_1", 600, 600);
91 hh_3d->Draw("iso");
92
93 // Draw the 2D projections of the 3D p.d.f.
94 TCanvas *c2 = new TCanvas("rf608_fitresultaspdf_2", "rf608_fitresultaspdf_2", 900, 600);
95 c2->Divide(3, 2);
96 c2->cd(1);
97 gPad->SetLeftMargin(0.15);
98 hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4);
99 hh_mean_sigmag2->Draw("surf3");
100 c2->cd(2);
101 gPad->SetLeftMargin(0.15);
102 hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4);
103 hh_sigmag2_frac->Draw("surf3");
104 c2->cd(3);
105 gPad->SetLeftMargin(0.15);
106 hh_mean_frac->GetZaxis()->SetTitleOffset(1.4);
107 hh_mean_frac->Draw("surf3");
108
109 // Draw the distributions of parameter points sampled from the p.d.f.
110 TH1 *tmp1 = d->createHistogram("mean,sigma_g2", 50, 50);
111 TH1 *tmp2 = d->createHistogram("sigma_g2,frac", 50, 50);
112 TH1 *tmp3 = d->createHistogram("mean,frac", 50, 50);
113
114 c2->cd(4);
115 gPad->SetLeftMargin(0.15);
116 tmp1->GetZaxis()->SetTitleOffset(1.4);
117 tmp1->Draw("lego3");
118 c2->cd(5);
119 gPad->SetLeftMargin(0.15);
120 tmp2->GetZaxis()->SetTitleOffset(1.4);
121 tmp2->Draw("lego3");
122 c2->cd(6);
123 gPad->SetLeftMargin(0.15);
124 tmp3->GetZaxis()->SetTitleOffset(1.4);
125 tmp3->Draw("lego3");
126}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define d(i)
Definition: RSha256.hxx:102
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:287
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition: RooAbsPdf.h:55
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:3342
TH1 * createHistogram(const char *varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const
Create and fill a ROOT histogram TH1, TH2 or TH3 with the values of this function for the variables w...
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:27
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
The 3-D histogram classes derived from the 1-D histogram classes.
Definition: TH3.h:31
RooCmdArg Save(Bool_t flag=kTRUE)
Double_t x[n]
Definition: legend1.C:17
return c2
Definition: legend2.C:14
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...