ROOT logo
//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #313
// 
// Working with parameterized ranges to define non-rectangular regions
// for fitting and integration
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

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


void rf313_paramranges()
{

  // C r e a t e   3 D   p d f 
  // -------------------------

  // Define observable (x,y,z)
  RooRealVar x("x","x",0,10) ;
  RooRealVar y("y","y",0,10) ;
  RooRealVar z("z","z",0,10) ;

  // Define 3 dimensional pdf
  RooRealVar z0("z0","z0",-0.1,1) ;
  RooPolynomial px("px","px",x,RooConst(0)) ;
  RooPolynomial py("py","py",y,RooConst(0)) ;
  RooPolynomial pz("pz","pz",z,z0) ;
  RooProdPdf pxyz("pxyz","pxyz",RooArgSet(px,py,pz)) ;



  // D e f i n e d   n o n - r e c t a n g u l a r   r e g i o n   R   i n   ( x , y , z ) 
  // -------------------------------------------------------------------------------------

  //
  // R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
  //

  // Construct range parameterized in "R" in y [ 0.1*x, 0.9*x ]
  RooFormulaVar ylo("ylo","0.1*x",x) ;
  RooFormulaVar yhi("yhi","0.9*x",x) ;
  y.setRange("R",ylo,yhi) ;

  // Construct parameterized ranged "R" in z [ 0, 0.1*y^2 ]
  RooFormulaVar zlo("zlo","0.0*y",y) ;
  RooFormulaVar zhi("zhi","0.1*y*y",y) ;
  z.setRange("R",zlo,zhi) ;



  // C a l c u l a t e   i n t e g r a l   o f   n o r m a l i z e d   p d f   i n   R 
  // ----------------------------------------------------------------------------------

  // Create integral over normalized pdf model over x,y,z in "R" region
  RooAbsReal* intPdf = pxyz.createIntegral(RooArgSet(x,y,z),RooArgSet(x,y,z),"R") ;

  // Plot value of integral as function of pdf parameter z0
  RooPlot* frame = z0.frame(Title("Integral of pxyz over x,y,z in region R")) ;
  intPdf->plotOn(frame) ;



  new TCanvas("rf313_paramranges","rf313_paramranges",600,600) ;
  gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;

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