ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #315
// 
// Marginizalization of multi-dimensional p.d.f.s through integration
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataHist.h"
#include "RooGaussian.h"
#include "RooProdPdf.h"
#include "RooPolyVar.h"
#include "TH1.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "RooNumIntConfig.h"
#include "RooConstVar.h"
using namespace RooFit ;


void rf315_projectpdf()
{
  // C r e a t e   p d f   m ( x , y )  =  g x ( x | y ) * g ( y )
  // --------------------------------------------------------------

  // Increase default precision of numeric integration
  // as this exercise has high sensitivity to numeric integration precision
  RooAbsPdf::defaultIntegratorConfig()->setEpsRel(1e-8) ;
  RooAbsPdf::defaultIntegratorConfig()->setEpsAbs(1e-8) ;

  // Create observables
  RooRealVar x("x","x",-5,5) ;
  RooRealVar y("y","y",-2,2) ;

  // Create function f(y) = a0 + a1*y
  RooRealVar a0("a0","a0",0) ;
  RooRealVar a1("a1","a1",-1.5,-3,1) ;
  RooPolyVar fy("fy","fy",y,RooArgSet(a0,a1)) ;

  // Create gaussx(x,f(y),sx)
  RooRealVar sigmax("sigmax","width of gaussian",0.5) ;
  RooGaussian gaussx("gaussx","Gaussian in x with shifting mean in y",x,fy,sigmax) ;  

  // Create gaussy(y,0,2)
  RooGaussian gaussy("gaussy","Gaussian in y",y,RooConst(0),RooConst(2)) ;

  // Create gaussx(x,sx|y) * gaussy(y)
  RooProdPdf model("model","gaussx(x|y)*gaussy(y)",gaussy,Conditional(gaussx,x)) ;
  


  // M a r g i n a l i z e   m ( x , y )   t o   m ( x ) 
  // ----------------------------------------------------

  // modelx(x) = Int model(x,y) dy
  RooAbsPdf* modelx = model.createProjection(y) ;



  // 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 .
  // ------------------------------------------------------------------------------------------

  // Sample 1000 events from modelx
  RooAbsData* data = modelx->generateBinned(x,1000) ;

  // Fit modelx to toy data
  modelx->fitTo(*data,Verbose()) ;

  // Plot modelx over data
  RooPlot* frame = x.frame(40) ;
  data->plotOn(frame) ;
  modelx->plotOn(frame) ;

  // Make 2D histogram of model(x,y)
  TH1* hh = model.createHistogram("x,y") ;
  hh->SetLineColor(kBlue) ;


  TCanvas* c = new TCanvas("rf315_projectpdf","rf315_projectpdf",800,400) ;
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.20) ; hh->GetZaxis()->SetTitleOffset(2.5) ; hh->Draw("surf") ;
  
}
 rf315_projectpdf.C:1
 rf315_projectpdf.C:2
 rf315_projectpdf.C:3
 rf315_projectpdf.C:4
 rf315_projectpdf.C:5
 rf315_projectpdf.C:6
 rf315_projectpdf.C:7
 rf315_projectpdf.C:8
 rf315_projectpdf.C:9
 rf315_projectpdf.C:10
 rf315_projectpdf.C:11
 rf315_projectpdf.C:12
 rf315_projectpdf.C:13
 rf315_projectpdf.C:14
 rf315_projectpdf.C:15
 rf315_projectpdf.C:16
 rf315_projectpdf.C:17
 rf315_projectpdf.C:18
 rf315_projectpdf.C:19
 rf315_projectpdf.C:20
 rf315_projectpdf.C:21
 rf315_projectpdf.C:22
 rf315_projectpdf.C:23
 rf315_projectpdf.C:24
 rf315_projectpdf.C:25
 rf315_projectpdf.C:26
 rf315_projectpdf.C:27
 rf315_projectpdf.C:28
 rf315_projectpdf.C:29
 rf315_projectpdf.C:30
 rf315_projectpdf.C:31
 rf315_projectpdf.C:32
 rf315_projectpdf.C:33
 rf315_projectpdf.C:34
 rf315_projectpdf.C:35
 rf315_projectpdf.C:36
 rf315_projectpdf.C:37
 rf315_projectpdf.C:38
 rf315_projectpdf.C:39
 rf315_projectpdf.C:40
 rf315_projectpdf.C:41
 rf315_projectpdf.C:42
 rf315_projectpdf.C:43
 rf315_projectpdf.C:44
 rf315_projectpdf.C:45
 rf315_projectpdf.C:46
 rf315_projectpdf.C:47
 rf315_projectpdf.C:48
 rf315_projectpdf.C:49
 rf315_projectpdf.C:50
 rf315_projectpdf.C:51
 rf315_projectpdf.C:52
 rf315_projectpdf.C:53
 rf315_projectpdf.C:54
 rf315_projectpdf.C:55
 rf315_projectpdf.C:56
 rf315_projectpdf.C:57
 rf315_projectpdf.C:58
 rf315_projectpdf.C:59
 rf315_projectpdf.C:60
 rf315_projectpdf.C:61
 rf315_projectpdf.C:62
 rf315_projectpdf.C:63
 rf315_projectpdf.C:64
 rf315_projectpdf.C:65
 rf315_projectpdf.C:66
 rf315_projectpdf.C:67
 rf315_projectpdf.C:68
 rf315_projectpdf.C:69
 rf315_projectpdf.C:70
 rf315_projectpdf.C:71
 rf315_projectpdf.C:72
 rf315_projectpdf.C:73
 rf315_projectpdf.C:74
 rf315_projectpdf.C:75
 rf315_projectpdf.C:76
 rf315_projectpdf.C:77
 rf315_projectpdf.C:78
 rf315_projectpdf.C:79
 rf315_projectpdf.C:80
 rf315_projectpdf.C:81
 rf315_projectpdf.C:82
 rf315_projectpdf.C:83
 rf315_projectpdf.C:84
 rf315_projectpdf.C:85
 rf315_projectpdf.C:86
 rf315_projectpdf.C:87
 rf315_projectpdf.C:88
 rf315_projectpdf.C:89
 rf315_projectpdf.C:90
 rf315_projectpdf.C:91
 rf315_projectpdf.C:92
 rf315_projectpdf.C:93
 rf315_projectpdf.C:94