Loading [MathJax]/extensions/tex2jax.js
Logo ROOT   6.16/01
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
rf308_normintegration2d.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #308
5///
6/// Examples on normalization of p.d.f.s,
7/// integration of p.d.fs, construction
8/// of cumulative distribution functions from p.d.f.s
9/// in two dimensions
10///
11/// \macro_image
12/// \macro_output
13/// \macro_code
14/// \author 07/2008 - Wouter Verkerke
15
16
17#include "RooRealVar.h"
18#include "RooGaussian.h"
19#include "RooConstVar.h"
20#include "RooProdPdf.h"
21#include "RooAbsReal.h"
22#include "RooPlot.h"
23#include "TCanvas.h"
24#include "TAxis.h"
25#include "TH1.h"
26using namespace RooFit ;
27
28
30{
31 // S e t u p m o d e l
32 // ---------------------
33
34 // Create observables x,y
35 RooRealVar x("x","x",-10,10) ;
36 RooRealVar y("y","y",-10,10) ;
37
38 // Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2)
39 RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ;
40 RooGaussian gy("gy","gy",y,RooConst(+2),RooConst(2)) ;
41
42 // Create gxy = gx(x)*gy(y)
43 RooProdPdf gxy("gxy","gxy",RooArgSet(gx,gy)) ;
44
45
46
47 // 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
48 // --------------------------------------------------------------------------------------------------
49
50 // Return 'raw' unnormalized value of gx
51 cout << "gxy = " << gxy.getVal() << endl ;
52
53 // Return value of gxy normalized over x _and_ y in range [-10,10]
54 RooArgSet nset_xy(x,y) ;
55 cout << "gx_Norm[x,y] = " << gxy.getVal(&nset_xy) << endl ;
56
57 // Create object representing integral over gx
58 // which is used to calculate gx_Norm[x,y] == gx / gx_Int[x,y]
59 RooAbsReal* igxy = gxy.createIntegral(RooArgSet(x,y)) ;
60 cout << "gx_Int[x,y] = " << igxy->getVal() << endl ;
61
62 // NB: it is also possible to do the following
63
64 // Return value of gxy normalized over x in range [-10,10] (i.e. treating y as parameter)
65 RooArgSet nset_x(x) ;
66 cout << "gx_Norm[x] = " << gxy.getVal(&nset_x) << endl ;
67
68 // Return value of gxy normalized over y in range [-10,10] (i.e. treating x as parameter)
69 RooArgSet nset_y(y) ;
70 cout << "gx_Norm[y] = " << gxy.getVal(&nset_y) << endl ;
71
72
73
74 // 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
75 // ----------------------------------------------------------------------------
76
77 // Define a range named "signal" in x from -5,5
78 x.setRange("signal",-5,5) ;
79 y.setRange("signal",-3,3) ;
80
81 // Create an integral of gxy_Norm[x,y] over x and y in range "signal"
82 // This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
83 // range named "signal"
84 RooAbsReal* igxy_sig = gxy.createIntegral(RooArgSet(x,y),NormSet(RooArgSet(x,y)),Range("signal")) ;
85 cout << "gx_Int[x,y|signal]_Norm[x,y] = " << igxy_sig->getVal() << endl ;
86
87
88
89
90 // 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
91 // -----------------------------------------------------------------------------------------------------
92
93 // Create the cumulative distribution function of gx
94 // i.e. calculate Int[-10,x] gx(x') dx'
95 RooAbsReal* gxy_cdf = gxy.createCdf(RooArgSet(x,y)) ;
96
97 // Plot cdf of gx versus x
98 TH1* hh_cdf = gxy_cdf->createHistogram("hh_cdf",x,Binning(40),YVar(y,Binning(40))) ;
99 hh_cdf->SetLineColor(kBlue) ;
100
101 new TCanvas("rf308_normintegration2d","rf308_normintegration2d",600,600) ;
102 gPad->SetLeftMargin(0.15) ; hh_cdf->GetZaxis()->SetTitleOffset(1.8) ;
103 hh_cdf->Draw("surf") ;
104
105}
@ kBlue
Definition: Rtypes.h:63
#define gPad
Definition: TVirtualPad.h:286
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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 *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
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:502
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 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 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
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg NormSet(const RooArgSet &nset)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooConstVar & RooConst(Double_t val)
Ta Range(0, 0, 1, 1)