rf315_projectpdf.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #315
00004 // 
00005 // Marginizalization of multi-dimensional p.d.f.s through integration
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataHist.h"
00018 #include "RooGaussian.h"
00019 #include "RooProdPdf.h"
00020 #include "RooPolyVar.h"
00021 #include "TH1.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 #include "RooNumIntConfig.h"
00026 #include "RooConstVar.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf315_projectpdf()
00031 {
00032   // C r e a t e   p d f   m ( x , y )  =  g x ( x | y ) * g ( y )
00033   // --------------------------------------------------------------
00034 
00035   // Increase default precision of numeric integration
00036   // as this exercise has high sensitivity to numeric integration precision
00037   RooAbsPdf::defaultIntegratorConfig()->setEpsRel(1e-8) ;
00038   RooAbsPdf::defaultIntegratorConfig()->setEpsAbs(1e-8) ;
00039 
00040   // Create observables
00041   RooRealVar x("x","x",-5,5) ;
00042   RooRealVar y("y","y",-2,2) ;
00043 
00044   // Create function f(y) = a0 + a1*y
00045   RooRealVar a0("a0","a0",0) ;
00046   RooRealVar a1("a1","a1",-1.5,-3,1) ;
00047   RooPolyVar fy("fy","fy",y,RooArgSet(a0,a1)) ;
00048 
00049   // Create gaussx(x,f(y),sx)
00050   RooRealVar sigmax("sigmax","width of gaussian",0.5) ;
00051   RooGaussian gaussx("gaussx","Gaussian in x with shifting mean in y",x,fy,sigmax) ;  
00052 
00053   // Create gaussy(y,0,2)
00054   RooGaussian gaussy("gaussy","Gaussian in y",y,RooConst(0),RooConst(2)) ;
00055 
00056   // Create gaussx(x,sx|y) * gaussy(y)
00057   RooProdPdf model("model","gaussx(x|y)*gaussy(y)",gaussy,Conditional(gaussx,x)) ;
00058   
00059 
00060 
00061   // M a r g i n a l i z e   m ( x , y )   t o   m ( x ) 
00062   // ----------------------------------------------------
00063 
00064   // modelx(x) = Int model(x,y) dy
00065   RooAbsPdf* modelx = model.createProjection(y) ;
00066 
00067 
00068 
00069   // U s e   m a r g i n a l i z e d   p . d . f .   a s   r e g u l a r   1 - D   p . d . f .
00070   // ------------------------------------------------------------------------------------------
00071 
00072   // Sample 1000 events from modelx
00073   RooAbsData* data = modelx->generateBinned(x,1000) ;
00074 
00075   // Fit modelx to toy data
00076   modelx->fitTo(*data,Verbose()) ;
00077 
00078   // Plot modelx over data
00079   RooPlot* frame = x.frame(40) ;
00080   data->plotOn(frame) ;
00081   modelx->plotOn(frame) ;
00082 
00083   // Make 2D histogram of model(x,y)
00084   TH1* hh = model.createHistogram("x,y") ;
00085   hh->SetLineColor(kBlue) ;
00086 
00087 
00088   TCanvas* c = new TCanvas("rf315_projectpdf","rf315_projectpdf",800,400) ;
00089   c->Divide(2) ;
00090   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00091   c->cd(2) ; gPad->SetLeftMargin(0.20) ; hh->GetZaxis()->SetTitleOffset(2.5) ; hh->Draw("surf") ;
00092   
00093 }

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