Logo ROOT   6.14/05
Reference Guide
rf402_datahandling.C File Reference

Detailed Description

View in nbviewer Open in SWAN 'DATA AND CATEGORIES' RooFit tutorial macro #402

Tools for manipulation of (un)binned datasets

pict1_rf402_datahandling.C.png
Processing /mnt/build/workspace/root-makedoc-v614/rootspi/rdoc/src/v6-14-00-patches/tutorials/roofit/rf402_datahandling.C...
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
DataStore d (d)
Contains 1000 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) y = 31.607 L(0 - 40) "y"
3) c = Plus(idx = 1)
"c"
1) 0x260f820 RooRealVar:: x = 9 L(-10 - 10) "x"
2) 0x23ac470 RooRealVar:: y = 31.607 L(0 - 40) "y"
3) 0x2679380 RooCategory:: c = Plus(idx = 1)
"c"
1) 0x260f820 RooRealVar:: x = 8 L(-10 - 10) "x"
2) 0x23ac470 RooRealVar:: y = 30 L(0 - 40) "y"
3) 0x2679380 RooCategory:: c = Minus(idx = -1)
"c"
>> d1 has only columns x,c
DataStore d (d)
Contains 1000 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) c = Plus(idx = 1)
"c"
>> d2 has only column y
DataStore d (d)
Contains 1000 entries
Observables:
1) y = 31.607 L(0 - 40) "y"
>> d3 has only the points with y>5.17
DataStore d (d)
Contains 973 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) y = 31.607 L(0 - 40) "y"
3) c = Plus(idx = 1)
"c"
>> d4 has only columns x,c for data points with y>5.17
DataStore d (d)
Contains 973 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) c = Plus(idx = 1)
"c"
>> merge d2(y) with d1(x,c) to form d1(x,c,y)
DataStore d (d)
Contains 1000 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) c = Plus(idx = 1)
"c"
3) y = 31.607 L(0 - 40) "y"
>> append data points of d3 to d1
DataStore d (d)
Contains 1973 entries
Observables:
1) x = 9 L(-10 - 10) "x"
2) c = Plus(idx = 1)
"c"
3) y = 31.607 L(0 - 40) "y"
>> construct dh (binned) from d(unbinned) but only take the x and y dimensions,
>> the category 'c' will be projected in the filling process
DataStore dh (binned version of d)
Contains 100 entries
Observables:
1) x = 9 L(-10 - 10) B(10) "x"
2) y = 31.607 L(0 - 40) B(10) "y"
Binned Dataset dh (binned version of d)
Contains 100 bins with a total weight of 1000
Observables: 1) x = 9 L(-10 - 10) B(10) "x"
2) y = 31.607 L(0 - 40) B(10) "y"
>> number of bins in dh : 100
>> sum of weights in dh : 1000
>> integral over histogram: 8000
>> retrieving the properties of the bin enclosing coordinate (x,y) = (0.3,20.5)
bin center:
1) 0x2385430 RooRealVar:: x = 1 L(-10 - 10) B(10) "x"
2) 0x23f4960 RooRealVar:: y = 22 L(0 - 40) B(10) "y"
weight = 76
>> Creating 1-dimensional projection on y of dh for bins with x>0
DataStore dh (binned version of d)
Contains 10 entries
Observables:
1) y = 38 L(0 - 40) B(10) "y"
Binned Dataset dh (binned version of d)
Contains 10 bins with a total weight of 500
Observables: 1) y = 38 L(0 - 40) B(10) "y"
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 500 will supercede previous event count of 1000 for normalization of PDF projections
>> Persisting d via ROOT I/O
TFile** rf402_datahandling.root
TFile* rf402_datahandling.root
KEY: RooDataSet d;1 d
KEY: TProcessID ProcessID0;1 3920539e-de7d-11e8-90cd-e2c88e80beef
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooDataHist.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooCategory.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "TFile.h"
using namespace RooFit ;
// WVE Add reduction by range
void rf402_datahandling()
{
// Binned (RooDataHist) and unbinned datasets (RooDataSet) share
// many properties and inherit from a common abstract base class
// (RooAbsData), that provides an interface for all operations
// that can be performed regardless of the data format
RooRealVar x("x","x",-10,10) ;
RooRealVar y("y","y", 0, 40) ;
RooCategory c("c","c") ;
c.defineType("Plus",+1) ;
c.defineType("Minus",-1) ;
// B a s i c O p e r a t i o n s o n u n b i n n e d d a t a s e t s
// --------------------------------------------------------------
// RooDataSet is an unbinned dataset (a collection of points in N-dimensional space)
RooDataSet d("d","d",RooArgSet(x,y,c)) ;
// Unlike RooAbsArgs (RooAbsPdf,RooFormulaVar,....) datasets are not attached to
// the variables they are constructed from. Instead they are attached to an internal
// clone of the supplied set of arguments
// Fill d with dummy values
Int_t i ;
for (i=0 ; i<1000 ; i++) {
x = i/50 - 10 ;
y = sqrt(1.0*i) ;
c.setLabel((i%2)?"Plus":"Minus") ;
// We must explicitly refer to x,y,c here to pass the values because
// d is not linked to them (as explained above)
d.add(RooArgSet(x,y,c)) ;
}
d.Print("v") ;
cout << endl ;
// The get() function returns a pointer to the internal copy of the RooArgSet(x,y,c)
// supplied in the constructor
const RooArgSet* row = d.get() ;
row->Print("v") ;
cout << endl ;
// Get with an argument loads a specific data point in row and returns
// a pointer to row argset. get() always returns the same pointer, unless
// an invalid row number is specified. In that case a null ptr is returned
d.get(900)->Print("v") ;
cout << endl ;
// R e d u c i n g , A p p e n d i n g a n d M e r g i n g
// -------------------------------------------------------------
// The reduce() function returns a new dataset which is a subset of the original
cout << endl << ">> d1 has only columns x,c" << endl ;
RooDataSet* d1 = (RooDataSet*) d.reduce(RooArgSet(x,c)) ;
d1->Print("v") ;
cout << endl << ">> d2 has only column y" << endl ;
d2->Print("v") ;
cout << endl << ">> d3 has only the points with y>5.17" << endl ;
RooDataSet* d3 = (RooDataSet*) d.reduce("y>5.17") ;
d3->Print("v") ;
cout << endl << ">> d4 has only columns x,c for data points with y>5.17" << endl ;
RooDataSet* d4 = (RooDataSet*) d.reduce(RooArgSet(x,c),"y>5.17") ;
d4->Print("v") ;
// The merge() function adds two data set column-wise
cout << endl << ">> merge d2(y) with d1(x,c) to form d1(x,c,y)" << endl ;
d1->merge(d2) ;
d1->Print("v") ;
// The append() function addes two datasets row-wise
cout << endl << ">> append data points of d3 to d1" << endl ;
d1->append(*d3) ;
d1->Print("v") ;
// O p e r a t i o n s o n b i n n e d d a t a s e t s
// ---------------------------------------------------------
// A binned dataset can be constructed empty, from an unbinned dataset, or
// from a ROOT native histogram (TH1,2,3)
cout << ">> construct dh (binned) from d(unbinned) but only take the x and y dimensions," << endl
<< ">> the category 'c' will be projected in the filling process" << endl ;
// The binning of real variables (like x,y) is done using their fit range
// 'get/setRange()' and number of specified fit bins 'get/setBins()'.
// Category dimensions of binned datasets get one bin per defined category state
x.setBins(10) ;
y.setBins(10) ;
RooDataHist dh("dh","binned version of d",RooArgSet(x,y),d) ;
dh.Print("v") ;
RooPlot* yframe = y.frame(Bins(10),Title("Operations on binned datasets")) ;
dh.plotOn(yframe) ; // plot projection of 2D binned data on y
// Examine the statistics of a binned dataset
cout << ">> number of bins in dh : " << dh.numEntries() << endl ;
cout << ">> sum of weights in dh : " << dh.sum(kFALSE) << endl ;
cout << ">> integral over histogram: " << dh.sum(kTRUE) << endl ; // accounts for bin volume
// Locate a bin from a set of coordinates and retrieve its properties
x = 0.3 ; y = 20.5 ;
cout << ">> retrieving the properties of the bin enclosing coordinate (x,y) = (0.3,20.5) " << endl ;
cout << " bin center:" << endl ;
dh.get(RooArgSet(x,y))->Print("v") ; // load bin center coordinates in internal buffer
cout << " weight = " << dh.weight() << endl ; // return weight of last loaded coordinates
// Reduce the 2-dimensional binned dataset to a 1-dimensional binned dataset
//
// All reduce() methods are interfaced in RooAbsData. All reduction techniques
// demonstrated on unbinned datasets can be applied to binned datasets as well.
cout << ">> Creating 1-dimensional projection on y of dh for bins with x>0" << endl ;
RooDataHist* dh2 = (RooDataHist*) dh.reduce(y,"x>0") ;
dh2->Print("v") ;
// Add dh2 to yframe and redraw
dh2->plotOn(yframe,LineColor(kRed),MarkerColor(kRed)) ;
// S a v i n g a n d l o a d i n g f r o m f i l e
// -------------------------------------------------------
// Datasets can be persisted with ROOT I/O
cout << endl << ">> Persisting d via ROOT I/O" << endl ;
TFile f("rf402_datahandling.root","RECREATE") ;
d.Write() ;
f.ls() ;
// To read back in future session:
// > TFile f("rf402_datahandling.root") ;
// > RooDataSet* d = (RooDataSet*) f.FindObject("d") ;
new TCanvas("rf402_datahandling","rf402_datahandling",600,600) ;
gPad->SetLeftMargin(0.15) ; yframe->GetYaxis()->SetTitleOffset(1.4) ; yframe->Draw() ;
}
Author
07/2008 - Wouter Verkerke

Definition in file rf402_datahandling.C.