rf308_normintegration2d.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #308
00004 // 
00005 // Examples on normalization of p.d.f.s,
00006 // integration of p.d.fs, construction
00007 // of cumulative distribution functions from p.d.f.s
00008 // in two dimensions
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 //
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooGaussian.h"
00019 #include "RooConstVar.h"
00020 #include "RooProdPdf.h"
00021 #include "RooAbsReal.h"
00022 #include "RooPlot.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "TH1.h"
00026 using namespace RooFit ;
00027 
00028 
00029 void rf308_normintegration2d()
00030 {
00031   // S e t u p   m o d e l 
00032   // ---------------------
00033 
00034   // Create observables x,y
00035   RooRealVar x("x","x",-10,10) ;
00036   RooRealVar y("y","y",-10,10) ;
00037 
00038   // Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2) 
00039   RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ;
00040   RooGaussian gy("gy","gy",y,RooConst(+2),RooConst(2)) ;
00041 
00042   // Create gxy = gx(x)*gy(y)
00043   RooProdPdf gxy("gxy","gxy",RooArgSet(gx,gy)) ;
00044 
00045 
00046 
00047   // 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
00048   // --------------------------------------------------------------------------------------------------
00049 
00050   // Return 'raw' unnormalized value of gx
00051   cout << "gxy = " << gxy.getVal() << endl ;
00052   
00053   // Return value of gxy normalized over x _and_ y in range [-10,10]
00054   RooArgSet nset_xy(x,y) ;
00055   cout << "gx_Norm[x,y] = " << gxy.getVal(&nset_xy) << endl ;
00056 
00057   // Create object representing integral over gx
00058   // which is used to calculate  gx_Norm[x,y] == gx / gx_Int[x,y]
00059   RooAbsReal* igxy = gxy.createIntegral(RooArgSet(x,y)) ;
00060   cout << "gx_Int[x,y] = " << igxy->getVal() << endl ;
00061 
00062   // NB: it is also possible to do the following
00063 
00064   // Return value of gxy normalized over x in range [-10,10] (i.e. treating y as parameter)
00065   RooArgSet nset_x(x) ;
00066   cout << "gx_Norm[x] = " << gxy.getVal(&nset_x) << endl ;
00067 
00068   // Return value of gxy normalized over y in range [-10,10] (i.e. treating x as parameter)
00069   RooArgSet nset_y(y) ;
00070   cout << "gx_Norm[y] = " << gxy.getVal(&nset_y) << endl ;
00071 
00072 
00073 
00074   // 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
00075   // ----------------------------------------------------------------------------
00076 
00077   // Define a range named "signal" in x from -5,5
00078   x.setRange("signal",-5,5) ;
00079   y.setRange("signal",-3,3) ;
00080   
00081   // Create an integral of gxy_Norm[x,y] over x and y in range "signal"
00082   // This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
00083   // range named "signal"
00084   RooAbsReal* igxy_sig = gxy.createIntegral(RooArgSet(x,y),NormSet(RooArgSet(x,y)),Range("signal")) ;
00085   cout << "gx_Int[x,y|signal]_Norm[x,y] = " << igxy_sig->getVal() << endl ;
00086 
00087 
00088 
00089 
00090   // 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
00091   // -----------------------------------------------------------------------------------------------------
00092 
00093   // Create the cumulative distribution function of gx
00094   // i.e. calculate Int[-10,x] gx(x') dx'
00095   RooAbsReal* gxy_cdf = gxy.createCdf(RooArgSet(x,y)) ;
00096   
00097   // Plot cdf of gx versus x
00098   TH1* hh_cdf = gxy_cdf->createHistogram("hh_cdf",x,Binning(40),YVar(y,Binning(40))) ;
00099   hh_cdf->SetLineColor(kBlue) ;
00100 
00101   new TCanvas("rf308_normintegration2d","rf308_normintegration2d",600,600) ;
00102   gPad->SetLeftMargin(0.15) ; hh_cdf->GetZaxis()->SetTitleOffset(1.8) ; 
00103   hh_cdf->Draw("surf") ;
00104 
00105 }

Generated on Tue Jul 5 15:45:05 2011 for ROOT_528-00b_version by  doxygen 1.5.1