00001 ///////////////////////////////////////////////////////////////////////// 00002 // 00003 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #110 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 one dimension 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 "RooAbsReal.h" 00021 #include "RooPlot.h" 00022 #include "TCanvas.h" 00023 #include "TAxis.h" 00024 using namespace RooFit ; 00025 00026 00027 void rf110_normintegration() 00028 { 00029 // S e t u p m o d e l 00030 // --------------------- 00031 00032 // Create observables x,y 00033 RooRealVar x("x","x",-10,10) ; 00034 00035 // Create p.d.f. gaussx(x,-2,3) 00036 RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ; 00037 00038 00039 // 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 00040 // -------------------------------------------------------------------------------------------------- 00041 00042 // Return 'raw' unnormalized value of gx 00043 cout << "gx = " << gx.getVal() << endl ; 00044 00045 // Return value of gx normalized over x in range [-10,10] 00046 RooArgSet nset(x) ; 00047 cout << "gx_Norm[x] = " << gx.getVal(&nset) << endl ; 00048 00049 // Create object representing integral over gx 00050 // which is used to calculate gx_Norm[x] == gx / gx_Int[x] 00051 RooAbsReal* igx = gx.createIntegral(x) ; 00052 cout << "gx_Int[x] = " << igx->getVal() << endl ; 00053 00054 00055 // 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 00056 // ---------------------------------------------------------------------------- 00057 00058 // Define a range named "signal" in x from -5,5 00059 x.setRange("signal",-5,5) ; 00060 00061 // Create an integral of gx_Norm[x] over x in range "signal" 00062 // This is the fraction of of p.d.f. gx_Norm[x] which is in the 00063 // range named "signal" 00064 RooAbsReal* igx_sig = gx.createIntegral(x,NormSet(x),Range("signal")) ; 00065 cout << "gx_Int[x|signal]_Norm[x] = " << igx_sig->getVal() << endl ; 00066 00067 00068 00069 // 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 00070 // ----------------------------------------------------------------------------------------------------- 00071 00072 // Create the cumulative distribution function of gx 00073 // i.e. calculate Int[-10,x] gx(x') dx' 00074 RooAbsReal* gx_cdf = gx.createCdf(x) ; 00075 00076 // Plot cdf of gx versus x 00077 RooPlot* frame = x.frame(Title("c.d.f of Gaussian p.d.f")) ; 00078 gx_cdf->plotOn(frame) ; 00079 00080 00081 00082 // Draw plot on canvas 00083 new TCanvas("rf110_normintegration","rf110_normintegration",600,600) ; 00084 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; 00085 frame->Draw() ; 00086 00087 00088 }