ROOT  6.06/09
Reference Guide
rf308_normintegration2d.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #308
4 //
5 // Examples on normalization of p.d.f.s,
6 // integration of p.d.fs, construction
7 // of cumulative distribution functions from p.d.f.s
8 // in two dimensions
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooGaussian.h"
19 #include "RooProdPdf.h"
20 #include "RooAbsReal.h"
21 #include "RooPlot.h"
22 #include "TCanvas.h"
23 #include "TH1.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic308 : public RooFitTestUnit
28 {
29 public:
30  TestBasic308(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Normalization of p.d.f.s in 2D",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // S e t u p m o d e l
34  // ---------------------
35 
36  // Create observables x,y
37  RooRealVar x("x","x",-10,10) ;
38  RooRealVar y("y","y",-10,10) ;
39 
40  // Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2)
41  RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ;
42  RooGaussian gy("gy","gy",y,RooConst(+2),RooConst(2)) ;
43 
44  // Create gxy = gx(x)*gy(y)
45  RooProdPdf gxy("gxy","gxy",RooArgSet(gx,gy)) ;
46 
47 
48 
49  // 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
50  // --------------------------------------------------------------------------------------------------
51 
52  // Return 'raw' unnormalized value of gx
53  regValue(gxy.getVal(),"rf308_gxy") ;
54 
55  // Return value of gxy normalized over x _and_ y in range [-10,10]
56  RooArgSet nset_xy(x,y) ;
57  regValue(gxy.getVal(&nset_xy),"rf308_gx_Norm[x,y]") ;
58 
59  // Create object representing integral over gx
60  // which is used to calculate gx_Norm[x,y] == gx / gx_Int[x,y]
61  RooAbsReal* igxy = gxy.createIntegral(RooArgSet(x,y)) ;
62  regValue(igxy->getVal(),"rf308_gx_Int[x,y]") ;
63 
64  // NB: it is also possible to do the following
65 
66  // Return value of gxy normalized over x in range [-10,10] (i.e. treating y as parameter)
67  RooArgSet nset_x(x) ;
68  regValue(gxy.getVal(&nset_x),"rf308_gx_Norm[x]") ;
69 
70  // Return value of gxy normalized over y in range [-10,10] (i.e. treating x as parameter)
71  RooArgSet nset_y(y) ;
72  regValue(gxy.getVal(&nset_y),"rf308_gx_Norm[y]") ;
73 
74 
75 
76  // 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
77  // ----------------------------------------------------------------------------
78 
79  // Define a range named "signal" in x from -5,5
80  x.setRange("signal",-5,5) ;
81  y.setRange("signal",-3,3) ;
82 
83  // Create an integral of gxy_Norm[x,y] over x and y in range "signal"
84  // This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
85  // range named "signal"
86  RooAbsReal* igxy_sig = gxy.createIntegral(RooArgSet(x,y),NormSet(RooArgSet(x,y)),Range("signal")) ;
87  regValue(igxy_sig->getVal(),"rf308_gx_Int[x,y|signal]_Norm[x,y]") ;
88 
89 
90 
91  // 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
92  // -----------------------------------------------------------------------------------------------------
93 
94  // Create the cumulative distribution function of gx
95  // i.e. calculate Int[-10,x] gx(x') dx'
96  RooAbsReal* gxy_cdf = gxy.createCdf(RooArgSet(x,y)) ;
97 
98  // Plot cdf of gx versus x
99  TH1* hh_cdf = gxy_cdf->createHistogram("hh_cdf",x,Binning(40),YVar(y,Binning(40))) ;
100 
101  regTH(hh_cdf,"rf308_cdf") ;
102 
103  delete igxy_sig ;
104  delete igxy ;
105  delete gxy_cdf ;
106 
107  return kTRUE ;
108  }
109 } ;
RooCmdArg NormSet(const RooArgSet &nset)
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Double_t x[n]
Definition: legend1.C:17
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
bool verbose
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
Double_t y[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:80
RooConstVar & RooConst(Double_t val)
const Bool_t kTRUE
Definition: Rtypes.h:91
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:503
RooCmdArg Binning(const RooAbsBinning &binning)