From $ROOTSYS/tutorials/roofit/rf108_plotbinning.C

/////////////////////////////////////////////////////////////////////////
//
// 'BASIC FUNCTIONALITY' RooFit tutorial macro #108
// 
// Plotting unbinned data with alternate and variable binnings
//
// 
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussModel.h"
#include "RooDecay.h"
#include "RooBMixDecay.h"
#include "RooCategory.h"
#include "RooBinning.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit ;


void rf108_plotbinning()
{

  // S e t u p   m o d e l 
  // ---------------------

  // Build a B decay p.d.f with mixing
  RooRealVar dt("dt","dt",-20,20) ;
  RooRealVar dm("dm","dm",0.472) ;
  RooRealVar tau("tau","tau",1.547) ;
  RooRealVar w("w","mistag rate",0.1) ;
  RooRealVar dw("dw","delta mistag rate",0.) ;

  RooCategory mixState("mixState","B0/B0bar mixing state") ;
  mixState.defineType("mixed",-1) ;
  mixState.defineType("unmixed",1) ;
  RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
  tagFlav.defineType("B0",1) ;
  tagFlav.defineType("B0bar",-1) ;

  // Build a gaussian resolution model
  RooRealVar dterr("dterr","dterr",0.1,1.0) ;
  RooRealVar bias1("bias1","bias1",0) ;
  RooRealVar sigma1("sigma1","sigma1",0.1) ;  
  RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;

  // Construct Bdecay (x) gauss
  RooBMixDecay bmix("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,RooBMixDecay::DoubleSided) ;


  // S a m p l e   d a t a   f r o m   m o d e l
  // --------------------------------------------

  // Sample 2000 events in (dt,mixState,tagFlav) from bmix
  RooDataSet *data = bmix.generate(RooArgSet(dt,mixState,tagFlav),2000) ;



  // S h o w   d t   d i s t r i b u t i o n   w i t h   c u s t o m   b i n n i n g
  // -------------------------------------------------------------------------------

  // Make plot of dt distribution of data in range (-15,15) with fine binning for dt>0 and coarse binning for dt<0

  // Create binning object with range (-15,15)
  RooBinning tbins(-15,15) ;

  // Add 60 bins with uniform spacing in range (-15,0)
  tbins.addUniform(60,-15,0) ;

  // Add 15 bins with uniform spacing in range (0,15)
  tbins.addUniform(15,0,15) ;
  
  // Make plot with specified binning
  RooPlot* dtframe = dt.frame(Range(-15,15),Title("dt distribution with custom binning")) ;
  data->plotOn(dtframe,Binning(tbins)) ;
  bmix.plotOn(dtframe) ;

  // NB: Note that bin density for each bin is adjusted to that of default frame binning as shown
  // in Y axis label (100 bins --> Events/0.4*Xaxis-dim) so that all bins represent a consistent density distribution
  
  
  // S h o w   m i x s t a t e   a s y m m e t r y  w i t h   c u s t o m   b i n n i n g
  // ------------------------------------------------------------------------------------

  // Make plot of dt distribution of data asymmetry in 'mixState' with variable binning 

  // Create binning object with range (-10,10)
  RooBinning abins(-10,10) ;

  // Add boundaries at 0, (-1,1), (-2,2), (-3,3), (-4,4) and (-6,6)
  abins.addBoundary(0) ;
  abins.addBoundaryPair(1) ;
  abins.addBoundaryPair(2) ;
  abins.addBoundaryPair(3) ;
  abins.addBoundaryPair(4) ;
  abins.addBoundaryPair(6) ;
  
  // Create plot frame in dt
  RooPlot* aframe = dt.frame(Range(-10,10),Title("mixState asymmetry distribution with custom binning")) ;

  // Plot mixState asymmetry of data with specified customg binning
  data->plotOn(aframe,Asymmetry(mixState),Binning(abins)) ;

  // Plot corresponding property of p.d.f
  bmix.plotOn(aframe,Asymmetry(mixState)) ;

  // Adjust vertical range of plot to sensible values for an asymmetry
  aframe->SetMinimum(-1.1) ;
  aframe->SetMaximum(1.1) ;

  // NB: For asymmetry distributions no density corrects are needed (and are thus not applied)


  // Draw plots on canvas
  TCanvas* c = new TCanvas("rf108_plotbinning","rf108_plotbinning",800,400) ;
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; dtframe->GetYaxis()->SetTitleOffset(1.6) ; dtframe->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; aframe->GetYaxis()->SetTitleOffset(1.6) ; aframe->Draw() ;
  
}
 rf108_plotbinning.C:1
 rf108_plotbinning.C:2
 rf108_plotbinning.C:3
 rf108_plotbinning.C:4
 rf108_plotbinning.C:5
 rf108_plotbinning.C:6
 rf108_plotbinning.C:7
 rf108_plotbinning.C:8
 rf108_plotbinning.C:9
 rf108_plotbinning.C:10
 rf108_plotbinning.C:11
 rf108_plotbinning.C:12
 rf108_plotbinning.C:13
 rf108_plotbinning.C:14
 rf108_plotbinning.C:15
 rf108_plotbinning.C:16
 rf108_plotbinning.C:17
 rf108_plotbinning.C:18
 rf108_plotbinning.C:19
 rf108_plotbinning.C:20
 rf108_plotbinning.C:21
 rf108_plotbinning.C:22
 rf108_plotbinning.C:23
 rf108_plotbinning.C:24
 rf108_plotbinning.C:25
 rf108_plotbinning.C:26
 rf108_plotbinning.C:27
 rf108_plotbinning.C:28
 rf108_plotbinning.C:29
 rf108_plotbinning.C:30
 rf108_plotbinning.C:31
 rf108_plotbinning.C:32
 rf108_plotbinning.C:33
 rf108_plotbinning.C:34
 rf108_plotbinning.C:35
 rf108_plotbinning.C:36
 rf108_plotbinning.C:37
 rf108_plotbinning.C:38
 rf108_plotbinning.C:39
 rf108_plotbinning.C:40
 rf108_plotbinning.C:41
 rf108_plotbinning.C:42
 rf108_plotbinning.C:43
 rf108_plotbinning.C:44
 rf108_plotbinning.C:45
 rf108_plotbinning.C:46
 rf108_plotbinning.C:47
 rf108_plotbinning.C:48
 rf108_plotbinning.C:49
 rf108_plotbinning.C:50
 rf108_plotbinning.C:51
 rf108_plotbinning.C:52
 rf108_plotbinning.C:53
 rf108_plotbinning.C:54
 rf108_plotbinning.C:55
 rf108_plotbinning.C:56
 rf108_plotbinning.C:57
 rf108_plotbinning.C:58
 rf108_plotbinning.C:59
 rf108_plotbinning.C:60
 rf108_plotbinning.C:61
 rf108_plotbinning.C:62
 rf108_plotbinning.C:63
 rf108_plotbinning.C:64
 rf108_plotbinning.C:65
 rf108_plotbinning.C:66
 rf108_plotbinning.C:67
 rf108_plotbinning.C:68
 rf108_plotbinning.C:69
 rf108_plotbinning.C:70
 rf108_plotbinning.C:71
 rf108_plotbinning.C:72
 rf108_plotbinning.C:73
 rf108_plotbinning.C:74
 rf108_plotbinning.C:75
 rf108_plotbinning.C:76
 rf108_plotbinning.C:77
 rf108_plotbinning.C:78
 rf108_plotbinning.C:79
 rf108_plotbinning.C:80
 rf108_plotbinning.C:81
 rf108_plotbinning.C:82
 rf108_plotbinning.C:83
 rf108_plotbinning.C:84
 rf108_plotbinning.C:85
 rf108_plotbinning.C:86
 rf108_plotbinning.C:87
 rf108_plotbinning.C:88
 rf108_plotbinning.C:89
 rf108_plotbinning.C:90
 rf108_plotbinning.C:91
 rf108_plotbinning.C:92
 rf108_plotbinning.C:93
 rf108_plotbinning.C:94
 rf108_plotbinning.C:95
 rf108_plotbinning.C:96
 rf108_plotbinning.C:97
 rf108_plotbinning.C:98
 rf108_plotbinning.C:99
 rf108_plotbinning.C:100
 rf108_plotbinning.C:101
 rf108_plotbinning.C:102
 rf108_plotbinning.C:103
 rf108_plotbinning.C:104
 rf108_plotbinning.C:105
 rf108_plotbinning.C:106
 rf108_plotbinning.C:107
 rf108_plotbinning.C:108
 rf108_plotbinning.C:109
 rf108_plotbinning.C:110
 rf108_plotbinning.C:111
 rf108_plotbinning.C:112
 rf108_plotbinning.C:113
 rf108_plotbinning.C:114
 rf108_plotbinning.C:115
 rf108_plotbinning.C:116
 rf108_plotbinning.C:117
 rf108_plotbinning.C:118
 rf108_plotbinning.C:119
 rf108_plotbinning.C:120
 rf108_plotbinning.C:121
 rf108_plotbinning.C:122
 rf108_plotbinning.C:123
 rf108_plotbinning.C:124
 rf108_plotbinning.C:125
 rf108_plotbinning.C:126
 rf108_plotbinning.C:127
 rf108_plotbinning.C:128
 rf108_plotbinning.C:129