rf308_normintegration2d.cxx

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 "RooProdPdf.h"
00020 #include "RooAbsReal.h"
00021 #include "RooPlot.h"
00022 #include "TCanvas.h"
00023 #include "TH1.h"
00024 using namespace RooFit ;
00025 
00026 
00027 class TestBasic308 : public RooFitTestUnit
00028 {
00029 public: 
00030   TestBasic308(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Normalization of p.d.f.s in 2D",refFile,writeRef,verbose) {} ;
00031   Bool_t testCode() {
00032 
00033   // S e t u p   m o d e l 
00034   // ---------------------
00035 
00036   // Create observables x,y
00037   RooRealVar x("x","x",-10,10) ;
00038   RooRealVar y("y","y",-10,10) ;
00039 
00040   // Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2) 
00041   RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ;
00042   RooGaussian gy("gy","gy",y,RooConst(+2),RooConst(2)) ;
00043 
00044   // Create gxy = gx(x)*gy(y)
00045   RooProdPdf gxy("gxy","gxy",RooArgSet(gx,gy)) ;
00046 
00047 
00048 
00049   // 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
00050   // --------------------------------------------------------------------------------------------------
00051 
00052   // Return 'raw' unnormalized value of gx
00053   regValue(gxy.getVal(),"rf308_gxy") ;
00054   
00055   // Return value of gxy normalized over x _and_ y in range [-10,10]
00056   RooArgSet nset_xy(x,y) ;
00057   regValue(gxy.getVal(&nset_xy),"rf308_gx_Norm[x,y]") ;
00058 
00059   // Create object representing integral over gx
00060   // which is used to calculate  gx_Norm[x,y] == gx / gx_Int[x,y]
00061   RooAbsReal* igxy = gxy.createIntegral(RooArgSet(x,y)) ;
00062   regValue(igxy->getVal(),"rf308_gx_Int[x,y]") ;
00063 
00064   // NB: it is also possible to do the following
00065 
00066   // Return value of gxy normalized over x in range [-10,10] (i.e. treating y as parameter)
00067   RooArgSet nset_x(x) ;
00068   regValue(gxy.getVal(&nset_x),"rf308_gx_Norm[x]") ;
00069 
00070   // Return value of gxy normalized over y in range [-10,10] (i.e. treating x as parameter)
00071   RooArgSet nset_y(y) ;
00072   regValue(gxy.getVal(&nset_y),"rf308_gx_Norm[y]") ;
00073 
00074 
00075 
00076   // 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
00077   // ----------------------------------------------------------------------------
00078 
00079   // Define a range named "signal" in x from -5,5
00080   x.setRange("signal",-5,5) ;
00081   y.setRange("signal",-3,3) ;
00082   
00083   // Create an integral of gxy_Norm[x,y] over x and y in range "signal"
00084   // This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
00085   // range named "signal"
00086   RooAbsReal* igxy_sig = gxy.createIntegral(RooArgSet(x,y),NormSet(RooArgSet(x,y)),Range("signal")) ;
00087   regValue(igxy_sig->getVal(),"rf308_gx_Int[x,y|signal]_Norm[x,y]") ;
00088 
00089 
00090 
00091   // 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
00092   // -----------------------------------------------------------------------------------------------------
00093 
00094   // Create the cumulative distribution function of gx
00095   // i.e. calculate Int[-10,x] gx(x') dx'
00096   RooAbsReal* gxy_cdf = gxy.createCdf(RooArgSet(x,y)) ;
00097   
00098   // Plot cdf of gx versus x
00099   TH1* hh_cdf = gxy_cdf->createHistogram("hh_cdf",x,Binning(40),YVar(y,Binning(40))) ;
00100   
00101   regTH(hh_cdf,"rf308_cdf") ;
00102 
00103   delete igxy_sig ;
00104   delete igxy ;
00105   delete gxy_cdf ;
00106 
00107   return kTRUE ;
00108   }
00109 } ;

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