From $ROOTSYS/tutorials/roofit/rf305_condcorrprod.C

/////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #305
// 
// Multi-dimensional p.d.f.s with conditional p.d.fs in product
// 
// pdf = gauss(x,f(y),sx | y ) * gauss(y,ms,sx)    with f(y) = a0 + a1*y
// 
//
// 07/2008 - Wouter Verkerke 
//
/////////////////////////////////////////////////////////////////////////

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



void rf305_condcorrprod()
{
  // C r e a t e   c o n d i t i o n a l   p d f   g x ( x | y ) 
  // -----------------------------------------------------------

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

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

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



  // C r e a t e   p d f   g y ( y ) 
  // -----------------------------------------------------------

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



  // C r e a t e   p r o d u c t   g x ( x | y ) * g y ( y )
  // -------------------------------------------------------

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



  // S a m p l e ,   f i t   a n d   p l o t   p r o d u c t   p d f
  // ---------------------------------------------------------------

  // Generate 1000 events in x and y from model
  RooDataSet *data = model.generate(RooArgSet(x,y),10000) ;

  // Plot x distribution of data and projection of model on x = Int(dy) model(x,y)
  RooPlot* xframe = x.frame() ;
  data->plotOn(xframe) ;
  model.plotOn(xframe) ; 

  // Plot x distribution of data and projection of model on y = Int(dx) model(x,y)
  RooPlot* yframe = y.frame() ;
  data->plotOn(yframe) ;
  model.plotOn(yframe) ; 

  // Make two-dimensional plot in x vs y
  TH1* hh_model = model.createHistogram("hh_model",x,Binning(50),YVar(y,Binning(50))) ;
  hh_model->SetLineColor(kBlue) ;



  // Make canvas and draw RooPlots
  TCanvas *c = new TCanvas("rf305_condcorrprod","rf05_condcorrprod",1200, 400);
  c->Divide(3);
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.6) ; xframe->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; yframe->GetYaxis()->SetTitleOffset(1.6) ; yframe->Draw() ;
  c->cd(3) ; gPad->SetLeftMargin(0.20) ; hh_model->GetZaxis()->SetTitleOffset(2.5) ; hh_model->Draw("surf") ;

}



 rf305_condcorrprod.C:1
 rf305_condcorrprod.C:2
 rf305_condcorrprod.C:3
 rf305_condcorrprod.C:4
 rf305_condcorrprod.C:5
 rf305_condcorrprod.C:6
 rf305_condcorrprod.C:7
 rf305_condcorrprod.C:8
 rf305_condcorrprod.C:9
 rf305_condcorrprod.C:10
 rf305_condcorrprod.C:11
 rf305_condcorrprod.C:12
 rf305_condcorrprod.C:13
 rf305_condcorrprod.C:14
 rf305_condcorrprod.C:15
 rf305_condcorrprod.C:16
 rf305_condcorrprod.C:17
 rf305_condcorrprod.C:18
 rf305_condcorrprod.C:19
 rf305_condcorrprod.C:20
 rf305_condcorrprod.C:21
 rf305_condcorrprod.C:22
 rf305_condcorrprod.C:23
 rf305_condcorrprod.C:24
 rf305_condcorrprod.C:25
 rf305_condcorrprod.C:26
 rf305_condcorrprod.C:27
 rf305_condcorrprod.C:28
 rf305_condcorrprod.C:29
 rf305_condcorrprod.C:30
 rf305_condcorrprod.C:31
 rf305_condcorrprod.C:32
 rf305_condcorrprod.C:33
 rf305_condcorrprod.C:34
 rf305_condcorrprod.C:35
 rf305_condcorrprod.C:36
 rf305_condcorrprod.C:37
 rf305_condcorrprod.C:38
 rf305_condcorrprod.C:39
 rf305_condcorrprod.C:40
 rf305_condcorrprod.C:41
 rf305_condcorrprod.C:42
 rf305_condcorrprod.C:43
 rf305_condcorrprod.C:44
 rf305_condcorrprod.C:45
 rf305_condcorrprod.C:46
 rf305_condcorrprod.C:47
 rf305_condcorrprod.C:48
 rf305_condcorrprod.C:49
 rf305_condcorrprod.C:50
 rf305_condcorrprod.C:51
 rf305_condcorrprod.C:52
 rf305_condcorrprod.C:53
 rf305_condcorrprod.C:54
 rf305_condcorrprod.C:55
 rf305_condcorrprod.C:56
 rf305_condcorrprod.C:57
 rf305_condcorrprod.C:58
 rf305_condcorrprod.C:59
 rf305_condcorrprod.C:60
 rf305_condcorrprod.C:61
 rf305_condcorrprod.C:62
 rf305_condcorrprod.C:63
 rf305_condcorrprod.C:64
 rf305_condcorrprod.C:65
 rf305_condcorrprod.C:66
 rf305_condcorrprod.C:67
 rf305_condcorrprod.C:68
 rf305_condcorrprod.C:69
 rf305_condcorrprod.C:70
 rf305_condcorrprod.C:71
 rf305_condcorrprod.C:72
 rf305_condcorrprod.C:73
 rf305_condcorrprod.C:74
 rf305_condcorrprod.C:75
 rf305_condcorrprod.C:76
 rf305_condcorrprod.C:77
 rf305_condcorrprod.C:78
 rf305_condcorrprod.C:79
 rf305_condcorrprod.C:80
 rf305_condcorrprod.C:81
 rf305_condcorrprod.C:82
 rf305_condcorrprod.C:83
 rf305_condcorrprod.C:84
 rf305_condcorrprod.C:85
 rf305_condcorrprod.C:86
 rf305_condcorrprod.C:87
 rf305_condcorrprod.C:88
 rf305_condcorrprod.C:89
 rf305_condcorrprod.C:90
 rf305_condcorrprod.C:91
 rf305_condcorrprod.C:92
 rf305_condcorrprod.C:93
 rf305_condcorrprod.C:94
 rf305_condcorrprod.C:95
 rf305_condcorrprod.C:96
 rf305_condcorrprod.C:97
 rf305_condcorrprod.C:98
 rf305_condcorrprod.C:99
 rf305_condcorrprod.C:100