Logo ROOT   6.16/01
Reference Guide
rf608_fitresultaspdf.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #608
5///
6/// Representing the parabolic approximation of the fit as
7/// a multi-variate Gaussian on the parameters of the fitted p.d.f.
8///
9/// \macro_image
10/// \macro_output
11/// \macro_code
12/// \author 07/2008 - Wouter Verkerke
13
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
32
34{
35 // C r e a t e m o d e l a n d d a t a s e t
36 // -----------------------------------------------
37
38 // Observable
39 RooRealVar x("x","x",-20,20) ;
40
41 // Model (intentional strong correlations)
42 RooRealVar mean("mean","mean of g1 and g2",0,-1,1) ;
43 RooRealVar sigma_g1("sigma_g1","width of g1",2) ;
44 RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
45
46 RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,5.0) ;
47 RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
48
49 RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
50 RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
51
52 // Generate 1000 events
53 RooDataSet* data = model.generate(x,1000) ;
54
55
56 // F i t m o d e l t o d a t a
57 // ----------------------------------
58
59 RooFitResult* r = model.fitTo(*data,Save()) ;
60
61
62 // 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
63 // ------------------------------------------------------------------------------------
64
65 RooAbsPdf* parabPdf = r->createHessePdf(RooArgSet(frac,mean,sigma_g2)) ;
66
67
68 // 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
69 // -----------------------------------------------------------------------------
70
71 // Generate 100K points in the parameter space, sampled from the MVGaussian p.d.f.
72 RooDataSet* d = parabPdf->generate(RooArgSet(mean,sigma_g2,frac),100000) ;
73
74
75 // Sample a 3-D histogram of the p.d.f. to be visualized as an error ellipsoid using the GLISO draw option
76 TH3* hh_3d = (TH3*) parabPdf->createHistogram("mean,sigma_g2,frac",25,25,25) ;
77 hh_3d->SetFillColor(kBlue) ;
78
79
80 // Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s
81 // The integrations corresponding to these projections are performed analytically
82 // by the MV Gaussian p.d.f.
83 RooAbsPdf* pdf_sigmag2_frac = parabPdf->createProjection(mean) ;
84 RooAbsPdf* pdf_mean_frac = parabPdf->createProjection(sigma_g2) ;
85 RooAbsPdf* pdf_mean_sigmag2 = parabPdf->createProjection(frac) ;
86
87
88 // Make 2D plots of the 3 two-dimensional p.d.f. projections
89 TH2* hh_sigmag2_frac = (TH2*) pdf_sigmag2_frac->createHistogram("sigma_g2,frac",50,50) ;
90 TH2* hh_mean_frac = (TH2*) pdf_mean_frac->createHistogram("mean,frac",50,50) ;
91 TH2* hh_mean_sigmag2 = (TH2*) pdf_mean_sigmag2->createHistogram("mean,sigma_g2",50,50) ;
92 hh_mean_frac->SetLineColor(kBlue) ;
93 hh_sigmag2_frac->SetLineColor(kBlue) ;
94 hh_mean_sigmag2->SetLineColor(kBlue) ;
95
96
97 // Draw the 'sigar'
98 new TCanvas("rf608_fitresultaspdf_1","rf608_fitresultaspdf_1",600,600) ;
99 hh_3d->Draw("iso") ;
100
101 // Draw the 2D projections of the 3D p.d.f.
102 TCanvas* c2 = new TCanvas("rf608_fitresultaspdf_2","rf608_fitresultaspdf_2",900,600) ;
103 c2->Divide(3,2) ;
104 c2->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_sigmag2->Draw("surf3") ;
105 c2->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_sigmag2_frac->Draw("surf3") ;
106 c2->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_mean_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_frac->Draw("surf3") ;
107
108 // Draw the distributions of parameter points sampled from the p.d.f.
109 TH1* tmp1 = d->createHistogram("mean,sigma_g2",50,50) ;
110 TH1* tmp2 = d->createHistogram("sigma_g2,frac",50,50) ;
111 TH1* tmp3 = d->createHistogram("mean,frac",50,50) ;
112
113 c2->cd(4) ; gPad->SetLeftMargin(0.15) ; tmp1->GetZaxis()->SetTitleOffset(1.4) ; tmp1->Draw("lego3") ;
114 c2->cd(5) ; gPad->SetLeftMargin(0.15) ; tmp2->GetZaxis()->SetTitleOffset(1.4) ; tmp2->Draw("lego3") ;
115 c2->cd(6) ; gPad->SetLeftMargin(0.15) ; tmp3->GetZaxis()->SetTitleOffset(1.4) ; tmp3->Draw("lego3") ;
116
117}
118
ROOT::R::TRInterface & r
Definition: Object.C:4
#define d(i)
Definition: RSha256.hxx:102
@ kBlue
Definition: Rtypes.h:63
#define gPad
Definition: TVirtualPad.h:286
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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:56
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:2914
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
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:31
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 fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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: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:31
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:2974
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
Double_t x[n]
Definition: legend1.C:17
return c2
Definition: legend2.C:14
RooCmdArg Save(Bool_t flag=kTRUE)