Logo ROOT  
Reference Guide
rf308_normintegration2d.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4///
5/// Multidimensional models: normalization and integration of p.d.fs, construction of
6/// cumulative distribution functions from p.d.f.s in two dimensions
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 "RooGaussian.h"
17#include "RooConstVar.h"
18#include "RooProdPdf.h"
19#include "RooAbsReal.h"
20#include "RooPlot.h"
21#include "TCanvas.h"
22#include "TAxis.h"
23#include "TH1.h"
24using namespace RooFit;
25
27{
28 // S e t u p m o d e l
29 // ---------------------
30
31 // Create observables x,y
32 RooRealVar x("x", "x", -10, 10);
33 RooRealVar y("y", "y", -10, 10);
34
35 // Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2)
36 RooGaussian gx("gx", "gx", x, RooConst(-2), RooConst(3));
37 RooGaussian gy("gy", "gy", y, RooConst(+2), RooConst(2));
38
39 // Create gxy = gx(x)*gy(y)
40 RooProdPdf gxy("gxy", "gxy", RooArgSet(gx, gy));
41
42 // R e t r i e v e r a w & n o r m a l i z e d v a l u e s o f R o o F i t p . d . f . s
43 // --------------------------------------------------------------------------------------------------
44
45 // Return 'raw' unnormalized value of gx
46 cout << "gxy = " << gxy.getVal() << endl;
47
48 // Return value of gxy normalized over x _and_ y in range [-10,10]
49 RooArgSet nset_xy(x, y);
50 cout << "gx_Norm[x,y] = " << gxy.getVal(&nset_xy) << endl;
51
52 // Create object representing integral over gx
53 // which is used to calculate gx_Norm[x,y] == gx / gx_Int[x,y]
54 RooAbsReal *igxy = gxy.createIntegral(RooArgSet(x, y));
55 cout << "gx_Int[x,y] = " << igxy->getVal() << endl;
56
57 // NB: it is also possible to do the following
58
59 // Return value of gxy normalized over x in range [-10,10] (i.e. treating y as parameter)
60 RooArgSet nset_x(x);
61 cout << "gx_Norm[x] = " << gxy.getVal(&nset_x) << endl;
62
63 // Return value of gxy normalized over y in range [-10,10] (i.e. treating x as parameter)
64 RooArgSet nset_y(y);
65 cout << "gx_Norm[y] = " << gxy.getVal(&nset_y) << endl;
66
67 // I n t e g r a t e n o r m a l i z e d p d f o v e r s u b r a n g e
68 // ----------------------------------------------------------------------------
69
70 // Define a range named "signal" in x from -5,5
71 x.setRange("signal", -5, 5);
72 y.setRange("signal", -3, 3);
73
74 // Create an integral of gxy_Norm[x,y] over x and y in range "signal"
75 // This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
76 // range named "signal"
77 RooAbsReal *igxy_sig = gxy.createIntegral(RooArgSet(x, y), NormSet(RooArgSet(x, y)), Range("signal"));
78 cout << "gx_Int[x,y|signal]_Norm[x,y] = " << igxy_sig->getVal() << endl;
79
80 // C o n s t r u c t c u m u l a t i v e d i s t r i b u t i o n f u n c t i o n f r o m p d f
81 // -----------------------------------------------------------------------------------------------------
82
83 // Create the cumulative distribution function of gx
84 // i.e. calculate Int[-10,x] gx(x') dx'
85 RooAbsReal *gxy_cdf = gxy.createCdf(RooArgSet(x, y));
86
87 // Plot cdf of gx versus x
88 TH1 *hh_cdf = gxy_cdf->createHistogram("hh_cdf", x, Binning(40), YVar(y, Binning(40)));
89 hh_cdf->SetLineColor(kBlue);
90
91 new TCanvas("rf308_normintegration2d", "rf308_normintegration2d", 600, 600);
92 gPad->SetLeftMargin(0.15);
93 hh_cdf->GetZaxis()->SetTitleOffset(1.8);
94 hh_cdf->Draw("surf");
95}
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:287
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
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...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:560
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
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 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
RooCmdArg NormSet(const RooArgSet &nset)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooConstVar & RooConst(Double_t val)
RooCmdArg Binning(const RooAbsBinning &binning)
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Ta Range(0, 0, 1, 1)