ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf608_fitresultaspdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #608
4 ///
5 /// Representing the parabolic approximation of the fit as
6 /// a multi-variate Gaussian on the parameters of the fitted p.d.f.
7 ///
8 ///
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooConstVar.h"
20 #include "RooAddPdf.h"
21 #include "RooChebychev.h"
22 #include "RooFitResult.h"
23 #include "TCanvas.h"
24 #include "TAxis.h"
25 #include "RooPlot.h"
26 #include "TFile.h"
27 #include "TStyle.h"
28 #include "TH2.h"
29 #include "TH3.h"
30 
31 using namespace RooFit ;
32 
33 
34 void rf608_fitresultaspdf()
35 {
36  // C r e a t e m o d e l a n d d a t a s e t
37  // -----------------------------------------------
38 
39  // Observable
40  RooRealVar x("x","x",-20,20) ;
41 
42  // Model (intentional strong correlations)
43  RooRealVar mean("mean","mean of g1 and g2",0,-1,1) ;
44  RooRealVar sigma_g1("sigma_g1","width of g1",2) ;
45  RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
46 
47  RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,5.0) ;
48  RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
49 
50  RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
51  RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
52 
53  // Generate 1000 events
54  RooDataSet* data = model.generate(x,1000) ;
55 
56 
57  // F i t m o d e l t o d a t a
58  // ----------------------------------
59 
60  RooFitResult* r = model.fitTo(*data,Save()) ;
61 
62 
63  // 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
64  // ------------------------------------------------------------------------------------
65 
66  RooAbsPdf* parabPdf = r->createHessePdf(RooArgSet(frac,mean,sigma_g2)) ;
67 
68 
69  // 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
70  // -----------------------------------------------------------------------------
71 
72  // Generate 100K points in the parameter space, sampled from the MVGaussian p.d.f.
73  RooDataSet* d = parabPdf->generate(RooArgSet(mean,sigma_g2,frac),100000) ;
74 
75 
76  // Sample a 3-D histogram of the p.d.f. to be visualized as an error ellipsoid using the GLISO draw option
77  TH3* hh_3d = (TH3*) parabPdf->createHistogram("mean,sigma_g2,frac",25,25,25) ;
78  hh_3d->SetFillColor(kBlue) ;
79 
80 
81  // Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s
82  // The integrations corresponding to these projections are performed analytically
83  // by the MV Gaussian p.d.f.
84  RooAbsPdf* pdf_sigmag2_frac = parabPdf->createProjection(mean) ;
85  RooAbsPdf* pdf_mean_frac = parabPdf->createProjection(sigma_g2) ;
86  RooAbsPdf* pdf_mean_sigmag2 = parabPdf->createProjection(frac) ;
87 
88 
89  // Make 2D plots of the 3 two-dimensional p.d.f. projections
90  TH2* hh_sigmag2_frac = (TH2*) pdf_sigmag2_frac->createHistogram("sigma_g2,frac",50,50) ;
91  TH2* hh_mean_frac = (TH2*) pdf_mean_frac->createHistogram("mean,frac",50,50) ;
92  TH2* hh_mean_sigmag2 = (TH2*) pdf_mean_sigmag2->createHistogram("mean,sigma_g2",50,50) ;
93  hh_mean_frac->SetLineColor(kBlue) ;
94  hh_sigmag2_frac->SetLineColor(kBlue) ;
95  hh_mean_sigmag2->SetLineColor(kBlue) ;
96 
97 
98  // Draw the 'sigar'
100  gStyle->SetPalette(1) ;
101  new TCanvas("rf608_fitresultaspdf_1","rf608_fitresultaspdf_1",600,600) ;
102  hh_3d->Draw("gliso") ;
103 
104  // Draw the 2D projections of the 3D p.d.f.
105  TCanvas* c2 = new TCanvas("rf608_fitresultaspdf_2","rf608_fitresultaspdf_2",900,600) ;
106  c2->Divide(3,2) ;
107  c2->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_sigmag2->Draw("surf3") ;
108  c2->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_sigmag2_frac->Draw("surf3") ;
109  c2->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_mean_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_frac->Draw("surf3") ;
110 
111  // Draw the distributions of parameter points sampled from the p.d.f.
112  TH1* tmp1 = d->createHistogram("mean,sigma_g2",50,50) ;
113  TH1* tmp2 = d->createHistogram("sigma_g2,frac",50,50) ;
114  TH1* tmp3 = d->createHistogram("mean,frac",50,50) ;
115 
116  c2->cd(4) ; gPad->SetLeftMargin(0.15) ; tmp1->GetZaxis()->SetTitleOffset(1.4) ; tmp1->Draw("lego3") ;
117  c2->cd(5) ; gPad->SetLeftMargin(0.15) ; tmp2->GetZaxis()->SetTitleOffset(1.4) ; tmp2->Draw("lego3") ;
118  c2->cd(6) ; gPad->SetLeftMargin(0.15) ; tmp3->GetZaxis()->SetTitleOffset(1.4) ; tmp3->Draw("lego3") ;
119 
120 }
121 
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:245
tuple g2
Definition: multifit.py:39
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
RooAbsPdf * createHessePdf(const RooArgSet &params) const
Return a p.d.f that represents the fit result as a multi-variate probability densisty function on the...
tuple g1
Definition: multifit.py:38
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
friend class RooArgSet
Definition: RooAbsArg.h:469
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:2989
The 3-D histogram classes derived from the 1-D histogram classes.
Definition: TH3.h:35
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
void SetCanvasPreferGL(Bool_t prefer=kTRUE)
Definition: TStyle.h:337
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
The TH1 histogram class.
Definition: TH1.h:80
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 wi...
TH2F * createHistogram(const RooAbsRealLValue &var1, const RooAbsRealLValue &var2, const char *cuts="", const char *name="hist") const
Create a TH2F histogram of the distribution of the specified variable using this dataset.
TAxis * GetZaxis()
Definition: TH1.h:321
RooCmdArg Save(Bool_t flag=kTRUE)
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())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:1702
#define gPad
Definition: TVirtualPad.h:288
Definition: Rtypes.h:61
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition: TStyle.cxx:1445
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:297