ROOT logo

From $ROOTSYS/tutorials/roofit/rf311_rangeplot.C

//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #310
// 
// Projecting p.d.f and data ranges in continuous observables
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

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


void rf311_rangeplot()
{

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

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

  // Create signal pdf gauss(x)*gauss(y)*gauss(z) 
  RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
  RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
  RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
  RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;

  // Create background pdf poly(x)*poly(y)*poly(z) 
  RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
  RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
  RooPolynomial pz("pz","pz",z) ;
  RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;

  // Create composite pdf sig+bkg
  RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
  RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;

  RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;



  // P r o j e c t   p d f   a n d   d a t a   o n   x
  // -------------------------------------------------

  // Make plain projection of data and pdf on x observable
  RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
  data->plotOn(frame) ;
  model.plotOn(frame) ;
  


  // P r o j e c t   p d f   a n d   d a t a   o n   x   i n   s i g n a l   r a n g e
  // ----------------------------------------------------------------------------------
  
  // Define signal region in y and z observables
  y.setRange("sigRegion",-1,1) ;
  z.setRange("sigRegion",-1,1) ;

  // Make plot frame
  RooPlot* frame2 = x.frame(Title("Same projection on X in signal range of (Y,Z)"),Bins(40)) ;

  // Plot subset of data in which all observables are inside "sigRegion"
  // For observables that do not have an explicit "sigRegion" range defined (e.g. observable)
  // an implicit definition is used that is identical to the full range (i.e. [-5,5] for x)
  data->plotOn(frame2,CutRange("sigRegion")) ;

  // Project model on x, integrating projected observables (y,z) only in "sigRegion"
  model.plotOn(frame2,ProjectionRange("sigRegion")) ;



  TCanvas* c = new TCanvas("rf311_rangeplot","rf310_rangeplot",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.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;


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