Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ROOT::RDataFrame Class Reference

ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees, CSV's and other data formats.

In addition, multi-threading and other low-level optimisations allow users to exploit all the resources available on their machines completely transparently.
Skip to the class reference or keep reading for the user guide.

In a nutshell:

ROOT::EnableImplicitMT(); // Tell ROOT you want to go parallel
ROOT::RDataFrame d("myTree", "file_*.root"); // Interface to TTree and TChain
auto myHisto = d.Histo1D("Branch_A"); // This books the (lazy) filling of a histogram
myHisto->Draw(); // Event loop is run here, upon first access to a result
#define d(i)
Definition RSha256.hxx:102
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition TROOT.cxx:525

Calculations are expressed in terms of a type-safe functional chain of actions and transformations, RDataFrame takes care of their execution. The implementation automatically puts in place several low level optimisations such as multi-thread parallelisation and caching.

DOI

For the impatient user

You can directly see RDataFrame in action in our tutorials, in C++ or Python.

Table of Contents

Cheat sheet

These are the operations which can be performed with RDataFrame.

Transformations

Transformations are a way to manipulate the data.

Transformation Description
Define Creates a new column in the dataset. Example usages include adding a column that contains the invariant mass of a particle, or a selection of elements of an array (e.g. only the pts of "good" muons).
DefineSlot Same as Define, but the user-defined function must take an extra unsigned int slot as its first parameter. slot will take a different value, 0 to nThreads - 1, for each thread of execution. This is meant as a helper in writing thread-safe Define transformation when using RDataFrame after ROOT::EnableImplicitMT(). DefineSlot works just as well with single-thread execution: in that case slot will always be 0.
DefineSlotEntry Same as DefineSlot, but the entry number is passed in addition to the slot number. This is meant as a helper in case some dependency on the entry number needs to be honoured.
Filter Filter rows based on user-defined conditions.
Range Filter rows based on entry number (single-thread only).

Actions

Actions aggregate data into a result. Each one is described in more detail in the reference guide.

In the following, whenever we say an action "returns" something, we always mean it returns a smart pointer to it. Actions only act on events that pass all preceding filters.

Lazy actions only trigger the event loop when one of the results is accessed for the first time, making it easy to produce many different results in one event loop. Instant actions trigger the event loop instantly.

Lazy action Description
Aggregate Execute a user-defined accumulation operation on the processed column values.
Book Book execution of a custom action using a user-defined helper object.
Cache Caches in contiguous memory columns' entries. Custom columns can be cached as well, filtered entries are not cached. Users can specify which columns to save (default is all).
Count Return the number of events processed. Useful e.g. to get a quick count of the number of events passing a Filter.
Display Provides a printable representation of the dataset contents. The method returns a RDisplay instance which can be queried to get a compressed tabular representation on the standard output or a complete representation as a string.
Fill Fill a user-defined object with the values of the specified columns, as if by calling Obj.Fill(col1, col2, ...). \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Graph](classROOT_1_1RDF_1_1RInterface.html#a804b466ebdbddef5c7e3400cc6b89301) \ilinebr </td> <td class="markdownTableBodyNone"> Fills a TGraph with the two columns provided. If Multithread is enabled, the order of the points may not be the one expected, it is therefore suggested to sort if before drawing. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [Histo{1D,2D,3D}](classROOT_1_1RDF_1_1RInterface.html#a247ca3aeb7ce5b95015b7fae72983055) \ilinebr </td> <td class="markdownTableBodyNone"> Fill a {one,two,three}-dimensional histogram with the processed column values. \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Max](classROOT_1_1RDF_1_1RInterface.html#a057179b1e77599466a0b02200d5cd8c3) \ilinebr </td> <td class="markdownTableBodyNone"> Return the maximum of processed column values. If the type of the column is inferred, the return type isdouble, the type of the column otherwise. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [Mean](classROOT_1_1RDF_1_1RInterface.html#ade6b020284f2f4fe9d3b09246b5f376a) \ilinebr </td> <td class="markdownTableBodyNone"> Return the mean of processed column values. \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Min](classROOT_1_1RDF_1_1RInterface.html#a7005702189e601972b6d19ecebcdc80c) \ilinebr </td> <td class="markdownTableBodyNone"> Return the minimum of processed column values. If the type of the column is inferred, the return type isdouble, the type of the column otherwise. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [Profile{1D,2D}](classROOT_1_1RDF_1_1RInterface.html#a8ef7dc16b0e9f7bc9cfbe2d9e5de0cef) \ilinebr </td> <td class="markdownTableBodyNone"> Fill a {one,two}-dimensional profile with the column values that passed all filters. \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Reduce](classROOT_1_1RDF_1_1RInterface.html#a118e723ae29834df8f2a992ded347354) \ilinebr </td> <td class="markdownTableBodyNone"> Reduce (e.g. sum, merge) entries using the function (lambda, functor...) passed as argument. The function must have signatureT(T,T)whereTis the type of the column. Return the final result of the reduction operation. An optional parameter allows initialization of the result object to non-default values. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [Report](classROOT_1_1RDF_1_1RInterface.html#a94f322531dcb25beb8f53a602e5d6332) \ilinebr </td> <td class="markdownTableBodyNone"> Obtains statistics on how many entries have been accepted and rejected by the filters. See the section on [named filters](#named-filters-and-cutflow-reports) for a more detailed explanation. The method returns a RCutFlowReport instance which can be queried programmatically to get information about the effects of the individual cuts. \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Stats](https://root.cern/doc/master/classROOT_1_1RDF_1_1RInterface.html#a9e8fafb75abfa4faed4da18dcde01568) \ilinebr </td> <td class="markdownTableBodyNone"> Return a TStatistic object filled with the input columns. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [StdDev](classROOT_1_1RDF_1_1RInterface.html#a482c4e4f81fe1e421c016f89cd281572) \ilinebr </td> <td class="markdownTableBodyNone"> Return the unbiased standard deviation of the processed column values. \ilinebr </td> </tr> <tr class="markdownTableRowOdd"> <td class="markdownTableBodyNone"> [Sum](classROOT_1_1RDF_1_1RInterface.html#a61d03407459120df6749af43ed506891) \ilinebr </td> <td class="markdownTableBodyNone"> Return the sum of the values in the column. If the type of the column is inferred, the return type isdouble, the type of the column otherwise. \ilinebr </td> </tr> <tr class="markdownTableRowEven"> <td class="markdownTableBodyNone"> [Take](classROOT_1_1RDF_1_1RInterface.html#a4fd694773a2931b6b07737ddcd1e73b4) \ilinebr </td> <td class="markdownTableBodyNone"> Extract a column from the dataset as a collection of values. If the type of the column is a C-style array, the type stored in the return container is aROOT::VecOps::RVec<T>` to guarantee the lifetime of the data involved.
Instant action Description
Foreach Execute a user-defined function on each entry. Users are responsible for the thread-safety of this lambda when executing with implicit multi-threading enabled.
ForeachSlot Same as Foreach, but the user-defined function must take an extra unsigned int slot as its first parameter. slot will take a different value, 0 to nThreads - 1, for each thread of execution. This is meant as a helper in writing thread-safe Foreach actions when using RDataFrame after ROOT::EnableImplicitMT(). ForeachSlot works just as well with single-thread execution: in that case slot will always be 0.
Snapshot Writes processed data-set to disk, in a new TTree and TFile. Custom columns can be saved as well, filtered entries are not saved. Users can specify which columns to save (default is all). Snapshot, by default, overwrites the output file if it already exists. Snapshot can be made lazy setting the appropriate flage in the snapshot options.

Other Operations

Operation Description
Alias Introduce an alias for a particular column name.
GetColumnNames Get the names of all the available columns of the dataset.
GetDefinedColumnNames Get the names of all the defined columns
GetColumnType Return the type of a given column as a string.
GetColumnTypeNamesList Return the list of type names of columns in the dataset.
GetFilterNames Get all the filters defined. If called on a root node, all filters will be returned. For any other node, only the filters upstream of that node.
SaveGraph Store the computation graph of an RDataFrame in graphviz format for easy inspection.
GetNRuns Get the number of event loops run by this RDataFrame instance.

Introduction

Users define their analysis as a sequence of operations to be performed on the data-frame object; the framework takes care of the management of the loop over entries as well as low-level details such as I/O and parallelisation. RDataFrame provides methods to perform most common operations required by ROOT analyses; at the same time, users can just as easily specify custom code that will be executed in the event loop.

RDataFrame is built with a modular and flexible workflow in mind, summarised as follows:

  1. build a data-frame object by specifying your data-set
  2. apply a series of transformations to your data
    1. filter (e.g. apply some cuts) or
    2. define a new column (e.g. the result of an expensive computation on columns)
  3. apply actions to the transformed data to produce results (e.g. fill a histogram)

Make sure to book all transformations and actions before you access the contents of any of the results: this lets RDataFrame accumulate work and then produce all results at the same time, upon first access to any of them.

The following table shows how analyses based on TTreeReader and TTree::Draw translate to RDataFrame. Follow the crash course to discover more idiomatic and flexible ways to express analyses with RDataFrame.

TTreeReader ROOT::RDataFrame
TTreeReader reader("myTree", file);
TTreeReaderValue<A_t> a(reader, "A");
TTreeReaderValue<B_t> b(reader, "B");
TTreeReaderValue<C_t> c(reader, "C");
while(reader.Next()) {
if(IsGoodEvent(*a, *b, *c))
DoStuff(*a, *b, *c);
}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
An interface for reading values stored in ROOT columnar datasets.
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition TTreeReader.h:44
Definition file.py:1
ROOT::RDataFrame d("myTree", file, {"A", "B", "C"});
d.Filter(IsGoodEvent).Foreach(DoStuff);
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, const ColumnNames_t &columns={}, std::string_view name="")
Append a filter to the call graph.
TTree::Draw ROOT::RDataFrame
auto t = file->Get<TTree>("myTree");
t->Draw("x", "y > 2");
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition TTree.h:428
ROOT::RDataFrame d("myTree", file);
auto h = d.Filter("y > 2").Histo1D("x");
#define h(i)
Definition RSha256.hxx:106

Crash course

All snippets of code presented in the crash course can be executed in the ROOT interpreter. Simply precede them with

using namespace ROOT; // RDataFrame's namespace
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...

which is omitted for brevity. The terms "column" and "branch" are used interchangeably.

Creating a RDataFrame

RDataFrame's constructor is where the user specifies the dataset and, optionally, a default set of columns that operations should work with. Here are the most common methods to construct a RDataFrame object:

// single file -- all ctors are equivalent
TFile *f = TFile::Open("file.root");
TTree *t = f.Get<TTree>("treeName");
RDataFrame d1("treeName", "file.root");
RDataFrame d2("treeName", f); // same as TTreeReader
RDataFrame d3(*t);
// multiple files -- all ctors are equivalent
TChain chain("myTree");
chain.Add("file1.root");
chain.Add("file2.root");
RDataFrame d4("myTree", {"file1.root", "file2.root"});
std::vector<std::string> files = {"file1.root", "file2.root"};
RDataFrame d5("myTree", files);
RDataFrame d6("myTree", "file*.root"); // the glob is passed as-is to TChain's constructor
RDataFrame d7(chain);
#define f(i)
Definition RSha256.hxx:104
A chain is a collection of files containing TTree objects.
Definition TChain.h:33
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:3997

Additionally, users can construct a RDataFrame with no data source by passing an integer number. This is the number of rows that will be generated by this RDataFrame.

RDataFrame d(10); // a RDF with 10 entries (and no columns/branches, for now)
d.Foreach([] { static int i = 0; std::cout << i++ << std::endl; }); // silly example usage: count to ten

This is useful to generate simple data-sets on the fly: the contents of each event can be specified with Define() (explained below). For example, we have used this method to generate Pythia events and write them to disk in parallel (with the Snapshot action).

For data sources other than TTrees and TChains, RDataFrame objects are constructed using ad-hoc factory functions (see e.g. MakeCsvDataFrame(), MakeSqliteDataFrame(), MakeArrowDataFrame()):

auto df = ROOT::RDF::MakeCsvDataFrame("input.csv");
// use df as usual
RDataFrame MakeCsvDataFrame(std::string_view fileName, bool readHeaders=true, char delimiter=',', Long64_t linesChunkSize=-1LL)
Factory method to create a CSV RDataFrame.
Definition RCsvDS.cxx:475

Filling a histogram

Let's now tackle a very common task, filling a histogram:

// Fill a TH1D with the "MET" branch
RDataFrame d("myTree", "file.root");
auto h = d.Histo1D("MET");
h->Draw();

The first line creates a RDataFrame associated to the TTree "myTree". This tree has a branch named "MET".

Histo1D() is an action; it returns a smart pointer (a RResultPtr, to be precise) to a TH1D histogram filled with the MET of all events. If the quantity stored in the branch is a collection (e.g. a vector or an array), the histogram is filled with all vector elements for each event.

You can use the objects returned by actions as if they were pointers to the desired results. There are many other possible actions, and all their results are wrapped in smart pointers; we'll see why in a minute.

Applying a filter

Let's say we want to cut over the value of branch "MET" and count how many events pass this cut. This is one way to do it:

RDataFrame d("myTree", "file.root");
auto c = d.Filter("MET > 4.").Count(); // computations booked, not run
std::cout << *c << std::endl; // computations run here, upon first access to the result

The filter string (which must contain a valid C++ expression) is applied to the specified branches for each event; the name and types of the columns are inferred automatically. The string expression is required to return a bool which signals whether the event passes the filter (true) or not (false).

You can think of your data as "flowing" through the chain of calls, being transformed, filtered and finally used to perform actions. Multiple Filter calls can be chained one after another.

Using string filters is nice for simple things, but they are limited to specifying the equivalent of a single return statement or the body of a lambda, so it's cumbersome to use strings with more complex filters. They also add a small runtime overhead, as ROOT needs to just-in-time compile the string into C++ code. When more freedom is required or runtime performance is very important, a C++ callable can be specified instead (a lambda in the following snippet, but it can be any kind of function or even a functor class), together with a list of branch names. This snippet is analogous to the one above:

RDataFrame d("myTree", "file.root");
auto metCut = [](double x) { return x > 4.; }; // a C++11 lambda function checking "x > 4"
auto c = d.Filter(metCut, {"MET"}).Count();
std::cout << *c << std::endl;
RResultPtr< ULong64_t > Count()
Return the number of entries processed (lazy action)
Double_t x[n]
Definition legend1.C:17

An example of a more complex filter expressed as a string containing C++ code is shown below

RDataFrame d("myTree", "file.root");
auto df = d.Define("p", "std::array<double, 4> p{px, py, pz}; return p;")
.Filter("double p2 = 0.0; for (auto&& x : p) p2 += x*x; return sqrt(p2) < 10.0;");

The code snippet above defines a column p that is a fixed-size array using the component column names and then filters on its magnitude by looping over its elements. It must be noted that the usage of strings to define columns like the one above is a major advantage when using PyROOT. However, only constants and data coming from other columns in the dataset can be involved in the code passed as a string. Local variables and functions cannot be used, since the interpreter will not know how to find them. When capturing local state is necessary, a C++ callable can be used.

More information on filters and how to use them to automatically generate cutflow reports can be found below.

Defining custom columns

Let's now consider the case in which "myTree" contains two quantities "x" and "y", but our analysis relies on a derived quantity z = sqrt(x*x + y*y). Using the Define transformation, we can create a new column in the data-set containing the variable "z":

RDataFrame d("myTree", "file.root");
auto sqrtSum = [](double x, double y) { return sqrt(x*x + y*y); };
auto zMean = d.Define("z", sqrtSum, {"x","y"}).Mean("z");
std::cout << *zMean << std::endl;
double sqrt(double)
RResultPtr< double > Mean(std::string_view columnName="")
Return the mean of processed column values (lazy action)
Double_t y[n]
Definition legend1.C:17

Define creates the variable "z" by applying sqrtSum to "x" and "y". Later in the chain of calls we refer to variables created with Define as if they were actual tree branches/columns, but they are evaluated on demand, at most once per event. As with filters, Define calls can be chained with other transformations to create multiple custom columns. Define and Filter transformations can be concatenated and intermixed at will.

As with filters, it is possible to specify new columns as string expressions. This snippet is analogous to the one above:

RDataFrame d("myTree", "file.root");
auto zMean = d.Define("z", "sqrt(x*x + y*y)").Mean("z");
std::cout << *zMean << std::endl;

Again the names of the branches used in the expression and their types are inferred automatically. The string must be valid C++ and is just-in-time compiled by the ROOT interpreter, cling – the process has a small runtime overhead.

Previously, when showing the different ways a RDataFrame can be created, we showed a constructor that only takes a number of entries a parameter. In the following example we show how to combine such an "empty" RDataFrame with Define transformations to create a data-set on the fly. We then save the generated data on disk using the Snapshot action.

RDataFrame d(100); // a RDF that will generate 100 entries (currently empty)
int x = -1;
auto d_with_columns = d.Define("x", [&x] { return ++x; })
.Define("xx", [&x] { return x*x; });
d_with_columns.Snapshot("myNewTree", "newfile.root");
RInterface< RDFDetail::RLoopManager, DS_t > Define(std::string_view name, F expression, const ColumnNames_t &columns={})
Creates a custom column.

This example is slightly more advanced than what we have seen so far: for starters, it makes use of lambda captures (a simple way to make external variables available inside the body of C++ lambdas) to act on the same variable x from both Define transformations. Secondly we have stored the transformed data-frame in a variable. This is always possible: at each point of the transformation chain, users can store the status of the data-frame for further use (more on this below).

You can read more about defining new columns here.

A graph composed of two branches, one starting with a filter and one with a define. The end point of a branch is always an action.

Running on a range of entries

It is sometimes necessary to limit the processing of the dataset to a range of entries. For this reason, the RDataFrame offers the concept of ranges as a node of the RDataFrame chain of transformations; this means that filters, columns and actions can be concatenated to and intermixed with Ranges. If a range is specified after a filter, the range will act exclusively on the entries passing the filter – it will not even count the other entries! The same goes for a Range hanging from another Range. Here are some commented examples:

RDataFrame d("myTree", "file.root");
// Here we store a data-frame that loops over only the first 30 entries in a variable
auto d30 = d.Range(30);
// This is how you pick all entries from 15 onwards
auto d15on = d.Range(15, 0);
// We can specify a stride too, in this case we pick an event every 3
auto d15each3 = d.Range(0, 15, 3);

Note that ranges are not available when multi-threading is enabled. More information on ranges is available here.

Executing multiple actions in the same event loop

As a final example let us apply two different cuts on branch "MET" and fill two different histograms with the "pt\_v" of the filtered events. By now, you should be able to easily understand what's happening:

RDataFrame d("treeName", "file.root");
auto h1 = d.Filter("MET > 10").Histo1D("pt_v");
auto h2 = d.Histo1D("pt_v");
h1->Draw(); // event loop is run once here
h2->Draw("SAME"); // no need to run the event loop again
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
TH1F * h1
Definition legend1.C:5

RDataFrame executes all above actions by running the event-loop only once. The trick is that actions are not executed at the moment they are called, but they are lazy, i.e. delayed until the moment one of their results is accessed through the smart pointer. At that time, the event loop is triggered and all results are produced simultaneously.

It is therefore good practice to declare all your transformations and actions before accessing their results, allowing RDataFrame to run the loop once and produce all results in one go.

Going parallel

Let's say we would like to run the previous examples in parallel on several cores, dividing events fairly between cores. The only modification required to the snippets would be the addition of this line before constructing the main data-frame object:

Simple as that. More details are given below.

Working with collections

RDataFrame reads collections as the special type RVec (e.g. a branch containing an array of floating point numbers can be read as a RVec<float>). C-style arrays (with variable or static size), std::vectors and most other collection types can be read this way. When reading ROOT data, column values of type RVec<T> perform no copy of the underlying array.

RVec is a container similar to std::vector but it offers a rich interface to operate on the array elements in a vectorised fashion, similar to Python's NumPy arrays.

For example, to fill a histogram with the pt of selected particles for each event, Define() can be used to create a column that contains the desired array elements as follows:

# h is filled with all the elements of `good_pts`, for each event
h = df.Define("good_pts", "pt[pt > 0]").Histo1D("good_pts")

Learn more on RVec.

Efficient analysis in Python

You can use RDataFrame in Python due to the dynamic C++/Python translation of PyROOT. In general, the interface is the same as for C++, a simple example follows.

df = ROOT.RDataFrame("myTree", "myFile.root")
sum = df.Filter("x > 10").Sum("y")
print(sum.GetValue())

Simple usage of efficient C++ code in Python

To perform more complex operations in the RDataFrame graph, e.g., in Filter and Define nodes, which don't fit into a simple expression string, you can just-in-time compile such functions directly in the Python script via the C++ interpreter cling. This approach has the advantage that you get the efficiency of compiled C++ code combined with the convenient workflow of a Python script. See the following snippet for an example of how to use a just-in-time-compiled C++ function from Python.

ROOT.gInterpreter.Declare("""
bool myFilter(float x) {
return x > 10;
}
""")
df = ROOT.RDataFrame("myTree", "myFile.root")
sum = df.Filter("myFilter(x)").Sum("y")
print(sum.GetValue())

To increase the performance even further, you can also pre-compile a C++ library with full code optimizations and load the function into the RDataFrame computation as follows.

ROOT.gSystem.Load("path/to/myLibrary.so") # Library with the myFilter function
ROOT.gInterpreter.Declare('#include "myLibrary.h"') # Header with the definition of the myFilter function
df = ROOT.RDataFrame("myTree", "myFile.root")
sum = df.Filter("myFilter(x)").Sum("y")
print(sum.GetValue())

Just-in-time compilation of Python callables with numba

ROOT also offers the option to compile Python callables with fundamental types and arrays thereof using numba and then using the function in RDataFrame from C++. The workflow requires the Python packages numba and cffi to be installed. See the following snippet for a simple example or the full tutorial here.

@ROOT.Numba.Declare(["float"], "bool")
def myFilter(x):
return x > 10
df = ROOT.RDataFrame("myTree", "myFile.root")
sum = df.Filter("Numba::myFilter(x)").Sum("y")
print(sum.GetValue())

Conversion to numpy arrays

Eventually, you probably would like to inspect the content of the RDataFrame or process the data further with functionality from Python libraries. For this purpose, we provide the AsNumpy function, which is able to provide you the columns of your RDataFrame as numpy arrays in Python. See a brief introduction below or a full tutorial here.

df = ROOT.RDataFrame("myTree", "myFile.root")
cols = df.Filter("x > 10").AsNumpy(["x", "y"])
print(cols["x"], cols["y"])

Distributed execution in Python

RDataFrame applications can be executed in parallel through distributed computing frameworks on a set of remote machines thanks to the Python package ROOT.RDF.Experimental.Distributed. This experimental, Python-only package allows to scale the optimized performance RDataFrame can achieve on a single machine to multiple nodes at the same time. It is designed so that different backends can be easily plugged in, currently supporting Apache Spark and soon also Dask. To make use of distributed RDataFrame, you only need to switch ROOT.RDataFrame with the backend-specific RDataFrame of your choice, for example:

import ROOT
# Point RDataFrame calls to the Spark specific RDataFrame
RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
# It still accepts the same constructor arguments as traditional RDataFrame
df = RDataFrame("mytree", "myfile.root")
# Continue the application with the traditional RDataFrame API
sum = df.Filter("x > 10").Sum("y")
h = df.Histo1D("x")
print(sum.GetValue())
h.Draw()

The main goal of this package is to support running any RDataFrame application distributedly. Nonetheless, not all RDataFrame operations currently work with this package. The subset that is currently available is:

  • AsNumpy
  • Count
  • Define
  • Fill
  • Filter
  • Graph
  • Histo[1,2,3]D
  • Max
  • Mean
  • Min
  • Profile[1,2,3]D
  • Snapshot
  • Sum

with support for more operations coming in the future. Data sources other than TTree and TChain (e.g. CSV, RNTuple) are currently not supported.

Connecting to a Spark cluster

In order to distribute the RDataFrame workload, you can connect to a Spark cluster you have access to through the official Spark API, then hook the connection instance to the distributed RDataFrame object like so:

import pyspark
import ROOT
# Create a SparkContext object with the right configuration for your Spark cluster
conf = SparkConf().setAppName(appName).setMaster(master)
sc = SparkContext(conf=conf)
# Point RDataFrame calls to the Spark specific RDataFrame
RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
# The Spark RDataFrame constructor accepts an optional "sparkcontext" parameter
# and it will distribute the application to the connected cluster
df = RDataFrame("mytree", "myfile.root", sparkcontext = sc)

If an instance of SparkContext is not provided, the default behaviour is to create one in the background for you.

Transformations

Filters

A filter is created through a call to Filter(f, columnList) or Filter(filterString). In the first overload, f can be a function, a lambda expression, a functor class, or any other callable object. It must return a bool signalling whether the event has passed the selection (true) or not (false). It should perform "read-only" operations on the columns, and should not have side-effects (e.g. modification of an external or static variable) to ensure correctness when implicit multi-threading is active. The second overload takes a string with a valid C++ expression in which column names are used as variable names (e.g. Filter("x[0] + x[1] > 0")). This is a convenience feature that comes with a certain runtime overhead: C++ code has to be generated on the fly from this expression before using it in the event loop. See the paragraph about "Just-in-time compilation" below for more information.

RDataFrame only evaluates filters when necessary: if multiple filters are chained one after another, they are executed in order and the first one returning false causes the event to be discarded and triggers the processing of the next entry. If multiple actions or transformations depend on the same filter, that filter is not executed multiple times for each entry: after the first access it simply serves a cached result.

Named filters and cutflow reports

An optional string parameter name can be passed to the Filter method to create a named filter. Named filters work as usual, but also keep track of how many entries they accept and reject.

Statistics are retrieved through a call to the Report method:

  • when Report is called on the main RDataFrame object, it returns a RResultPtr<RCutFlowReport> relative to all named filters declared up to that point
  • when called on a specific node (e.g. the result of a Define or Filter), it returns a RResultPtr<RCutFlowReport> relative all named filters in the section of the chain between the main RDataFrame and that node (included).

Stats are stored in the same order as named filters have been added to the graph, and refer to the latest event-loop that has been run using the relevant RDataFrame.

Ranges

When RDataFrame is not being used in a multi-thread environment (i.e. no call to EnableImplicitMT was made), Range transformations are available. These act very much like filters but instead of basing their decision on a filter expression, they rely on begin,end and stride parameters.

  • begin: initial entry number considered for this range.
  • end: final entry number (excluded) considered for this range. 0 means that the range goes until the end of the dataset.
  • stride: process one entry of the [begin, end) range every stride entries. Must be strictly greater than 0.

The actual number of entries processed downstream of a Range node will be (end - begin)/stride (or less if less entries than that are available).

Note that ranges act "locally", not based on the global entry count: Range(10,50) means "skip the first 10 entries that reach this node*, let the next 40 entries pass, then stop processing". If a range node hangs from a filter node, and the range has a begin parameter of 10, that means the range will skip the first 10 entries that pass the preceding filter.

Ranges allow "early quitting": if all branches of execution of a functional graph reached their end value of processed entries, the event-loop is immediately interrupted. This is useful for debugging and quick data explorations.

Custom columns

Custom columns are created by invoking Define(name, f, columnList). As usual, f can be any callable object (function, lambda expression, functor class...); it takes the values of the columns listed in columnList (a list of strings) as parameters, in the same order as they are listed in columnList. f must return the value that will be assigned to the temporary column.

A new variable is created called name, accessible as if it was contained in the dataset from subsequent transformations/actions.

Use cases include:

  • caching the results of complex calculations for easy and efficient multiple access
  • extraction of quantities of interest from complex objects
  • branch aliasing, i.e. changing the name of a branch

An exception is thrown if the name of the new column/branch is already in use for another branch in the TTree.

It is also possible to specify the quantity to be stored in the new temporary column as a C++ expression with the method Define(name, expression). For example this invocation

df.Define("pt", "sqrt(px*px + py*py)");

will create a new column called "pt" the value of which is calculated starting from the columns px and py. The system builds a just-in-time compiled function starting from the expression after having deduced the list of necessary branches from the names of the variables specified by the user.

Custom columns as function of slot and entry number

It is possible to create custom columns also as a function of the processing slot and entry numbers. The methods that can be invoked are:

Actions

Instant and lazy actions

Actions can be instant or lazy. Instant actions are executed as soon as they are called, while lazy actions are executed whenever the object they return is accessed for the first time. As a rule of thumb, actions with a return value are lazy, the others are instant.

Performance tips and parallel execution

As pointed out before in this document, RDataFrame can transparently perform multi-threaded event loops to speed up the execution of its actions. Users have to call ROOT::EnableImplicitMT() before constructing the RDataFrame object to indicate that it should take advantage of a pool of worker threads. Each worker thread processes a distinct subset of entries, and their partial results are merged before returning the final values to the user. There are no guarantees on the order in which threads will process the batches of entries. In particular, note that this means that, for multi-thread event loops, there is no guarantee on the order in which Snapshot will write entries: they could be scrambled with respect to the input dataset.

Warning
By default, RDataFrame will use as many threads as the hardware supports, using up all the resources on a machine. This might be undesirable on shared computing resources such as a batch cluster. Therefore, when running on shared computing resources, use replacing i with the number of CPUs/slots that were allocated for this job.

Thread-safety of user-defined expressions

RDataFrame operations such as Histo1D or Snapshot are guaranteed to work correctly in multi-thread event loops. User-defined expressions, such as strings or lambdas passed to Filter, Define, Foreach, Reduce or Aggregate will have to be thread-safe, i.e. it should be possible to call them concurrently from different threads.

Note that simple Filter and Define transformations will inherently satisfy this requirement: Filter/Define expressions will often be pure in the functional programming sense (no side-effects, no dependency on external state), which eliminates all risks of race conditions.

In order to facilitate writing of thread-safe operations, some RDataFrame features such as Foreach, Define or OnPartialResult offer thread-aware counterparts (ForeachSlot, DefineSlot, OnPartialResultSlot): their only difference is that they will pass an extra slot argument (an unsigned integer) to the user-defined expression. When calling user-defined code concurrently, RDataFrame guarantees that different threads will employ different values of the slot parameter, where slot will be a number between 0 and ROOT::GetThreadPoolSize() - 1. In other words, within a slot, computation runs sequentially and events are processed sequentially. Note that the same slot might be associated to different threads over the course of a single event loop, but two threads will never receive the same slot at the same time. This extra parameter might facilitate writing safe parallel code by having each thread write/modify a different processing slot*, e.g. a different element of a list. See here for an example usage of ForeachSlot.

Parallel execution of multiple RDataFrame event loops

A complex analysis may require multiple separate RDatFrame computation graphs to produce all desired results. This poses the challenge that the event loops of each computation graph can be parallelized, but the different loops run sequentially, one after the other. On many-core architectures it might be desirable to run different event loops concurrently to improve resource usage. ROOT::RDF::RunGraphs() allows running multiple RDataFrame event loops concurrently:

ROOT::RDataFrame df1("tree1", "f1.root");
ROOT::RDataFrame df2("tree2", "f2.root");
auto histo1 = df1.Histo1D("x");
auto histo2 = df2.Histo1D("y");
// just accessing result pointers, the event loops of separate RDataFrames run one after the other
histo1->Draw(); // runs first multi-thread event loop
histo2->Draw(); // runs second multi-thread event loop
// with ROOT::RDF::RunGraphs, event loops for separate computation graphs can run concurrently
ROOT::RDF::RunGraphs({histo1, histo2});
void RunGraphs(std::vector< RResultHandle > handles)
Trigger the event loop of multiple RDataFrames concurrently.

Performance profiling of RDataFrame applications

To obtain the maximum performance out of RDataFrame, make sure to avoid just-in-time compiled versions of transformations and actions if at all possible. For instance, Filter("x > 0") requires just-in-time compilation of the corresponding C++ logic, while the equivalent Filter([] { return x > 0; }, {"x"}) does not. Similarly, Histo1D("x") requires just-in-time compilation after the type of x is retrieved from the dataset, while Histo1D<float>("x") does not; the latter spelling should be preferred for performance-critical applications.

Python applications cannot easily specify template parameters or pass C++ callables to RDataFrame. See Efficient analysis in Python for possible ways to speed up hot paths in this case.

Just-in-time compilation happens once, right before starting an event loop. To reduce the runtime cost of this step, make sure to book all operations for all RDataFrame computation graphs before the first event loop is triggered: just-in-time compilation will happen once for all code required to be generated up to that point, also across different computation graphs.

Also make sure not to count the just-in-time compilation time (which happens once before the event loop and does not depend on the size of the dataset) as part of the event loop runtime (which scales with the size of the dataset). RDataFrame has an experimental logging feature that simplifies measuring the time spent in just-in-time compilation and in the event loop (as well as providing some more interesting information). It is activated like follows:

auto verbosity = ROOT::Experimental::RLogScopedVerbosity(ROOT::Detail::RDF::RDFLogChannel(), ROOT::Experimental::ELogLevel.kInfo);
const Int_t kInfo
Definition TError.h:46
ELogLevel
Kinds of diagnostics.
Definition RLogger.hxx:36

Memory usage

There are two reasons why RDataFrame may consume more memory than expected. Firstly, each result is duplicated for each worker thread, which e.g. in case of many (possibly multi-dimensional) histograms with fine binning can result in visible memory consumption during the event loop. The thread-local copies of the results are destroyed when the final result is produced. Secondly, just-in-time compilation of string expressions or non-templated actions (see the previous paragraph) causes Cling, ROOT's C++ interpreter, to allocate some memory for the generated code that is only released at the end of the application. This commonly results in memory usage creep in long-running applications that create many RDataFrames one after the other. Possible mitigations include creating and running each RDataFrame event loop in a sub-process, or booking all operations for all different RDataFrame computation graphs before the first event loop is triggered, so that the interpreter is invoked only once for all computation graphs.

More features

Here is a list of the most important features that have been omitted in the "Crash course" for brevity. You don't need to read all these to start using RDataFrame, but they are useful to save typing time and runtime.

Programmatically get the list of column names

The GetColumnsNames() method returns the list of valid column names for the dataset:

RDataFrame d("myTree", "file.root");
std::vector<std::string> colNames = d.GetColumnNames();

Callbacks

It's possible to schedule execution of arbitrary functions (callbacks) during the event loop. Callbacks can be used e.g. to inspect partial results of the analysis while the event loop is running, drawing a partially-filled histogram every time a certain number of new entries is processed, or event displaying a progress bar while the event loop runs.

For example one can draw an up-to-date version of a result histogram every 100 entries like this:

auto h = tdf.Histo1D("x");
TCanvas c("c","x hist");
h.OnPartialResult(100, [&c](TH1D &h_) { c.cd(); h_.Draw(); c.Update(); });
h->Draw(); // event loop runs here, this `Draw` is executed after the event loop is finished
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:618

Callbacks are registered to a RResultPtr and must be callables that takes a reference to the result type as argument and return nothing. RDataFrame will invoke registered callbacks passing partial action results as arguments to them (e.g. a histogram filled with a part of the selected events).

Read more on RResultPtr::OnPartialResult().

Default branch lists

When constructing a RDataFrame object, it is possible to specify a default column list for your analysis, in the usual form of a list of strings representing branch/column names. The default column list will be used as a fallback whenever a list specific to the transformation/action is not present. RDataFrame will take as many of these columns as needed, ignoring trailing extra names if present.

// use "b1" and "b2" as default branches
RDataFrame d1("myTree", "file.root", {"b1","b2"});
auto h = d1.Filter([](int b1, int b2) { return b1 > b2; }) // will act on "b1" and "b2"
.Histo1D(); // will act on "b1"
// just one default branch this time
RDataFrame d2("myTree", "file.root", {"b1"});
auto min = d2.Filter([](double b2) { return b2 > 0; }, {"b2"}) // we can still specify non-default branch lists
.Min(); // returns the minimum value of "b1" for the filtered entries
RResultPtr< RDFDetail::MinReturnType_t< T > > Min(std::string_view columnName="")
Return the minimum of processed column values (lazy action)
RResultPtr<::TH1D > Histo1D(const TH1DModel &model={"", "", 128u, 0., 0.}, std::string_view vName="")
Fill and return a one-dimensional histogram with the values of a column (lazy action)

Implicit Columns

Every instance of RDataFrame is created with two special columns called rdfentry_ and rdfslot_. The rdfentry_ column is an unsigned 64-bit integer holding the current entry number while rdfslot_ is an unsigned 32-bit integer holding the index of the current data processing slot. For backwards compatibility reasons, the names tdfentry_ and tdfslot_ are also accepted. These columns are not considered by operations such as Cache or Snapshot. The cached or snapshot data frame provides "its own" values for these columns which do not necessarily correspond to the ones of the mother data frame. This is most notably the case where filters are used before deriving a cached/persistified dataframe.

Note that in multi-thread event loops the values of rdfentry_ do not correspond to what would be the entry numbers of a TChain constructed over the same set of ROOT files, as the entries are processed in an unspecified order.

Just-in-time compilation: branch type inference and explicit declaration of branch types

C++ is a statically typed language: all types must be known at compile-time. This includes the types of the TTree branches we want to work on. For filters, temporary columns and some of the actions, branch types are deduced from the signature of the relevant filter function/temporary column expression/action function:

// here b1 is deduced to be `int` and b2 to be `double`
dataFrame.Filter([](int x, double y) { return x > 0 && y < 0.; }, {"b1", "b2"});

If we specify an incorrect type for one of the branches, an exception with an informative message will be thrown at runtime, when the branch value is actually read from the TTree: RDataFrame detects type mismatches. The same would happen if we swapped the order of "b1" and "b2" in the branch list passed to Filter.

Certain actions, on the other hand, do not take a function as argument (e.g. Histo1D), so we cannot deduce the type of the branch at compile-time. In this case **RDataFrame infers the type of the branch** from the TTree itself. This is why we never needed to specify the branch types for all actions in the above snippets.

When the branch type is not a common one such as int, double, char or float it is nonetheless good practice to specify it as a template parameter to the action itself, like this:

dataFrame.Histo1D("b1"); // OK, the type of "b1" is deduced at runtime
dataFrame.Min<MyNumber_t>("myObject"); // OK, "myObject" is deduced to be of type `MyNumber_t`

Deducing types at runtime requires the just-in-time compilation of the relevant actions, which has a small runtime overhead, so specifying the type of the columns as template parameters to the action is good practice when performance is a goal.

When deducing types at runtime, fundamental types are read as constant values, i.e. it is not possible to write to column values from Filters or Defines. This is typically perfectly fine and avoids certain common mistakes such as typing x = 0 rather than x == 0. Classes and other complex types are read by non-constant references to avoid copies and to permit calls to non-const member functions. Note that calling non-const member functions will often not be thread-safe.

Generic actions

RDataFrame strives to offer a comprehensive set of standard actions that can be performed on each event. At the same time, it allows users to execute arbitrary code (i.e. a generic action) inside the event loop through the Foreach and ForeachSlot actions.

Foreach(f, columnList) takes a function f (lambda expression, free function, functor...) and a list of columns, and executes f on those columns for each event. The function passed must return nothing (i.e. void). It can be used to perform actions that are not already available in the interface. For example, the following snippet evaluates the root mean square of column "b":

// Single-thread evaluation of RMS of column "b" using Foreach
double sumSq = 0.;
unsigned int n = 0;
RDataFrame d("bTree", bFilePtr);
d.Foreach([&sumSq, &n](double b) { ++n; sumSq += b*b; }, {"b"});
std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
const Int_t n
Definition legend1.C:16

When executing on multiple threads, users are responsible for the thread-safety of the expression passed to Foreach: each thread will execute the expression multiple times (once per entry) in an unspecified order. The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of rms; but this is probably too much head-scratch for such a simple operation.

ForeachSlot can help in this situation. It is an alternative version of Foreach for which the function takes an additional parameter besides the columns it should be applied to: an unsigned int slot parameter, where slot is a number indicating which thread (0, 1, 2 , ..., poolSize - 1) the function is being run in. More specifically, RDataFrame guarantees that ForeachSlot will invoke the user expression with different slot parameters for different concurrent executions (there is no guarantee that a certain slot number will always correspond to a given thread id, though). We can take advantage of ForeachSlot to evaluate a thread-safe root mean square of branch "b":

// Thread-safe evaluation of RMS of branch "b" using ForeachSlot
const unsigned int nSlots = ROOT::GetThreadPoolSize();
std::vector<double> sumSqs(nSlots, 0.);
std::vector<unsigned int> ns(nSlots, 0);
RDataFrame d("bTree", bFilePtr);
d.ForeachSlot([&sumSqs, &ns](unsigned int slot, double b) { sumSqs[slot] += b*b; ns[slot] += 1; }, {"b"});
double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
UInt_t GetThreadPoolSize()
Returns the size of ROOT's thread pool.
Definition TROOT.cxx:563

You see how we created one double variable for each thread in the pool, and later merged their results via std::accumulate.

Friend trees

Friend trees are supported by RDataFrame. In order to deal with friend trees with RDataFrame, the user is required to build the tree and its friends and instantiate a RDataFrame with it.

TTree t([...]);
TTree ft([...]);
t.AddFriend(ft, "myFriend");
auto f = d.Filter("myFriend.MyCol == 42");
virtual TFriendElement * AddFriend(const char *treename, const char *filename="")
Add a TFriendElement to the list of friends.
Definition TTree.cxx:1321

Reading file formats different from ROOT's

RDataFrame can be interfaced with RDataSources. The RDataSource interface defines an API that RDataFrame can use to read arbitrary data formats.

A concrete RDataSource implementation (i.e. a class that inherits from RDataSource and implements all of its pure methods) provides an adaptor that RDataFrame can leverage to read any kind of tabular data formats. RDataFrame calls into RDataSource to retrieve information about the data, retrieve (thread-local) readers or "cursors" for selected columns and to advance the readers to the desired data entry. Some predefined RDataSources are natively provided by ROOT such as the RCsvDS which allows to read comma separated files:

auto tdf = ROOT::RDF::MakeCsvDataFrame("MuRun2010B.csv");
auto filteredEvents =
tdf.Filter("Q1 * Q2 == -1")
.Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
auto h = filteredEvents.Histo1D("m");
h->Draw();

Call graphs (storing and reusing sets of transformations)

Sets of transformations can be stored as variables** and reused multiple times to create call graphs in which several paths of filtering/creation of columns are executed simultaneously; we often refer to this as "storing the state of the chain".

This feature can be used, for example, to create a temporary column once and use it in several subsequent filters or actions, or to apply a strict filter to the data-set before executing several other transformations and actions, effectively reducing the amount of events processed.

Let's try to make this clearer with a commented example:

// build the data-frame and specify a default column list
RDataFrame d(treeName, filePtr, {"var1", "var2", "var3"});
// apply a cut and save the state of the chain
auto filtered = d.Filter(myBigCut);
// plot branch "var1" at this point of the chain
auto h1 = filtered.Histo1D("var1");
// create a new branch "vec" with a vector extracted from a complex object (only for filtered entries)
// and save the state of the chain
auto newBranchFiltered = filtered.Define("vec", [](const Obj& o) { return o.getVector(); }, {"obj"});
// apply a cut and fill a histogram with "vec"
auto h2 = newBranchFiltered.Filter(cut1).Histo1D("vec");
// apply a different cut and fill a new histogram
auto h3 = newBranchFiltered.Filter(cut2).Histo1D("vec");
// Inspect results
h2->Draw(); // first access to an action result: run event-loop!
h3->Draw("SAME"); // event loop does not need to be run again here..
std::cout << "Entries in h1: " << h1->GetEntries() << std::endl; // ..or here
virtual Double_t GetEntries() const
Return the current number of entries.
Definition TH1.cxx:4386

RDataFrame detects when several actions use the same filter or the same temporary column, and only evaluates each filter or temporary column once per event, regardless of how many times that result is used down the call graph. Objects read from each column are built once and never copied, for maximum efficiency. When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated, so it might be advisable to put the strictest filters first in the chain.

Printing the computation graph

It is possible to print the computation graph from any node to obtain a dot representation either on the standard output or in a file.

Invoking the function ROOT::RDF::SaveGraph() on any node that is not the head node, the computation graph of the branch the node belongs to is printed. By using the head node, the entire computation graph is printed.

Following there is an example of usage:

// First, a sample computational graph is built
ROOT::RDataFrame df("tree", "f.root");
auto df2 = df.Define("x", []() { return 1; })
.Filter("col0 % 1 == col0")
.Filter([](int b1) { return b1 <2; }, {"cut1"})
.Define("y", []() { return 1; });
auto count = df2.Count();
// Prints the graph to the rd1.dot file in the current directory
ROOT::RDF::SaveGraph(rd1, "./mydot.dot");
// Prints the graph to standard output
std::string SaveGraph(NodeType node)
Create a graphviz representation of the dataframe computation graph, return it as a string.

RDataFrame variables as function arguments and return values

RDataFrame variables/nodes are relatively cheap to copy and it's possible to both pass them to (or move them into) functions and to return them from functions. However, in general each dataframe node will have a different C++ type, which includes all available compile-time information about what that node does. One way to cope with this complication is to use template functions and/or C++14 auto return types:

template <typename RDF>
auto ApplySomeFilters(RDF df)
{
return df.Filter("x > 0").Filter([](int y) { return y < 0; }, {"y"});
}

A possibly simpler, C++11-compatible alternative is to take advantage of the fact that any dataframe node can be converted to the common type ROOT::RDF::RNode:

// a function that conditionally adds a Range to a RDataFrame node.
RNode MaybeAddRange(RNode df, bool mustAddRange)
{
return mustAddRange ? df.Range(1) : df;
}
// use as :
auto maybeRangedDF = MaybeAddRange(df, true);

The conversion to ROOT::RDF::RNode is cheap, but it will introduce an extra virtual call during the RDataFrame event loop (in most cases, the resulting performance impact should be negligible).

As a final note, remember that RDataFrame actions do not return another dataframe, but a RResultPtr<T>, where T is the type of the result of the action.

Read more on this topic here.

Definition at line 42 of file RDataFrame.hxx.

Public Types

using ColumnNames_t = RDFDetail::ColumnNames_t
 

Public Member Functions

 RDataFrame (std::string_view treeName, ::TDirectory *dirPtr, const ColumnNames_t &defaultBranches={})
 
 RDataFrame (std::string_view treename, const std::vector< std::string > &filenames, const ColumnNames_t &defaultBranches={})
 Build the dataframe.
 
 RDataFrame (std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultBranches={})
 Build the dataframe.
 
 RDataFrame (std::unique_ptr< ROOT::RDF::RDataSource >, const ColumnNames_t &defaultBranches={})
 Build dataframe associated to datasource.
 
 RDataFrame (TTree &tree, const ColumnNames_t &defaultBranches={})
 Build the dataframe.
 
 RDataFrame (ULong64_t numEntries)
 Build a dataframe that generates numEntries entries.
 
- Public Member Functions inherited from ROOT::RDF::RInterface< RDFDetail::RLoopManager >
 RInterface (const RInterface &)=default
 Copy-ctor for RInterface.
 
 RInterface (const std::shared_ptr< RDFDetail::RLoopManager > &proxied)
 Only enabled when building a RInterface<RLoopManager>
 
 RInterface (RInterface &&)=default
 Move-ctor for RInterface.
 
RResultPtr< U > Aggregate (AccFun aggregator, MergeFun merger, std::string_view columnName, const U &aggIdentity)
 Execute a user-defined accumulation operation on the processed column values in each processing slot.
 
RResultPtr< U > Aggregate (AccFun aggregator, MergeFun merger, std::string_view columnName="")
 Execute a user-defined accumulation operation on the processed column values in each processing slot.
 
RInterface< RDFDetail::RLoopManager, DS_tAlias (std::string_view alias, std::string_view columnName)
 Allow to refer to a column with a different name.
 
RResultPtr< typename Helper::Result_t > Book (Helper &&helper, const ColumnNames_t &columns={})
 Book execution of a custom action using a user-defined helper object.
 
RInterface< RLoopManagerCache (const ColumnNames_t &columnList)
 Save selected columns in memory.
 
RInterface< RLoopManagerCache (const ColumnNames_t &columnList)
 Save selected columns in memory.
 
RInterface< RLoopManagerCache (std::initializer_list< std::string > columnList)
 Save selected columns in memory.
 
RInterface< RLoopManagerCache (std::string_view columnNameRegexp="")
 Save selected columns in memory.
 
RResultPtr< ULong64_tCount ()
 Return the number of entries processed (lazy action)
 
RInterface< RDFDetail::RLoopManager, DS_tDefine (std::string_view name, F expression, const ColumnNames_t &columns={})
 Creates a custom column.
 
RInterface< RDFDetail::RLoopManager, DS_tDefine (std::string_view name, std::string_view expression)
 Creates a custom column.
 
RInterface< RDFDetail::RLoopManager, DS_tDefineSlot (std::string_view name, F expression, const ColumnNames_t &columns={})
 Creates a custom column with a value dependent on the processing slot.
 
RInterface< RDFDetail::RLoopManager, DS_tDefineSlotEntry (std::string_view name, F expression, const ColumnNames_t &columns={})
 Creates a custom column with a value dependent on the processing slot and the current entry.
 
RResultPtr< RDisplayDisplay (const ColumnNames_t &columnList, const int &nRows=5)
 Provides a representation of the columns in the dataset.
 
RResultPtr< RDisplayDisplay (const ColumnNames_t &columnList, const int &nRows=5)
 Provides a representation of the columns in the dataset.
 
RResultPtr< RDisplayDisplay (std::initializer_list< std::string > columnList, const int &nRows=5)
 Provides a representation of the columns in the dataset.
 
RResultPtr< RDisplayDisplay (std::string_view columnNameRegexp="", const int &nRows=5)
 Provides a representation of the columns in the dataset.
 
RResultPtr< T > Fill (T &&model, const ColumnNames_t &columnList)
 Return an object of type T on which T::Fill will be called once per event (lazy action)
 
RResultPtr< T > Fill (T &&model, const ColumnNames_t &columnList)
 Return an object of type T on which T::Fill will be called once per event (lazy action)
 
RInterface< RDFDetail::RFilter< F, RDFDetail::RLoopManager >, DS_tFilter (F f, const ColumnNames_t &columns={}, std::string_view name="")
 Append a filter to the call graph.
 
RInterface< RDFDetail::RFilter< F, RDFDetail::RLoopManager >, DS_tFilter (F f, const std::initializer_list< std::string > &columns)
 Append a filter to the call graph.
 
RInterface< RDFDetail::RFilter< F, RDFDetail::RLoopManager >, DS_tFilter (F f, std::string_view name)
 Append a filter to the call graph.
 
RInterface< RDFDetail::RJittedFilter, DS_tFilter (std::string_view expression, std::string_view name="")
 Append a filter to the call graph.
 
void Foreach (F f, const ColumnNames_t &columns={})
 Execute a user-defined function on each entry (instant action)
 
void ForeachSlot (F f, const ColumnNames_t &columns={})
 Execute a user-defined function requiring a processing slot index on each entry (instant action)
 
ColumnNames_t GetColumnNames ()
 Returns the names of the available columns.
 
std::string GetColumnType (std::string_view column)
 Return the type of a given column as a string.
 
ColumnNames_t GetDefinedColumnNames ()
 Returns the names of the defined columns.
 
std::vector< std::string > GetFilterNames ()
 Returns the names of the filters created.
 
unsigned int GetNRuns () const
 Gets the number of event loops run.
 
unsigned int GetNSlots () const
 Gets the number of data processing slots.
 
RResultPtr<::TGraphGraph (std::string_view v1Name="", std::string_view v2Name="")
 Fill and return a graph (lazy action)
 
bool HasColumn (std::string_view columnName)
 Checks if a column is present in the dataset.
 
RResultPtr<::TH1DHisto1D (const TH1DModel &model, std::string_view vName, std::string_view wName)
 Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
 
RResultPtr<::TH1DHisto1D (const TH1DModel &model={"", "", 128u, 0., 0.})
 Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
 
RResultPtr<::TH1DHisto1D (const TH1DModel &model={"", "", 128u, 0., 0.}, std::string_view vName="")
 Fill and return a one-dimensional histogram with the values of a column (lazy action)
 
RResultPtr<::TH1DHisto1D (std::string_view vName)
 Fill and return a one-dimensional histogram with the values of a column (lazy action)
 
RResultPtr<::TH1DHisto1D (std::string_view vName, std::string_view wName)
 Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
 
RResultPtr<::TH2DHisto2D (const TH2DModel &model)
 
RResultPtr<::TH2DHisto2D (const TH2DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view wName)
 Fill and return a weighted two-dimensional histogram (lazy action)
 
RResultPtr<::TH2DHisto2D (const TH2DModel &model, std::string_view v1Name="", std::string_view v2Name="")
 Fill and return a two-dimensional histogram (lazy action)
 
RResultPtr<::TH3DHisto3D (const TH3DModel &model)
 
RResultPtr<::TH3DHisto3D (const TH3DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view v3Name, std::string_view wName)
 Fill and return a three-dimensional histogram (lazy action)
 
RResultPtr<::TH3DHisto3D (const TH3DModel &model, std::string_view v1Name="", std::string_view v2Name="", std::string_view v3Name="")
 Fill and return a three-dimensional histogram (lazy action)
 
RResultPtr< RDFDetail::MaxReturnType_t< T > > Max (std::string_view columnName="")
 Return the maximum of processed column values (lazy action)
 
RResultPtr< doubleMean (std::string_view columnName="")
 Return the mean of processed column values (lazy action)
 
RResultPtr< RDFDetail::MinReturnType_t< T > > Min (std::string_view columnName="")
 Return the minimum of processed column values (lazy action)
 
 operator RNode () const
 Cast any RDataFrame node to a common type ROOT::RDF::RNode.
 
RInterfaceoperator= (const RInterface &)=default
 Copy-assignment operator for RInterface.
 
RResultPtr<::TProfileProfile1D (const TProfile1DModel &model)
 
RResultPtr<::TProfileProfile1D (const TProfile1DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view wName)
 Fill and return a one-dimensional profile (lazy action)
 
RResultPtr<::TProfileProfile1D (const TProfile1DModel &model, std::string_view v1Name="", std::string_view v2Name="")
 Fill and return a one-dimensional profile (lazy action)
 
RResultPtr<::TProfile2DProfile2D (const TProfile2DModel &model)
 
RResultPtr<::TProfile2DProfile2D (const TProfile2DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view v3Name, std::string_view wName)
 Fill and return a two-dimensional profile (lazy action)
 
RResultPtr<::TProfile2DProfile2D (const TProfile2DModel &model, std::string_view v1Name="", std::string_view v2Name="", std::string_view v3Name="")
 Fill and return a two-dimensional profile (lazy action)
 
RInterface< RDFDetail::RRange< RDFDetail::RLoopManager >, DS_tRange (unsigned int begin, unsigned int end, unsigned int stride=1)
 Creates a node that filters entries based on range: [begin, end)
 
RInterface< RDFDetail::RRange< RDFDetail::RLoopManager >, DS_tRange (unsigned int end)
 Creates a node that filters entries based on range.
 
RResultPtr< T > Reduce (F f, std::string_view columnName, const T &redIdentity)
 Execute a user-defined reduce operation on the values of a column.
 
RResultPtr< T > Reduce (F f, std::string_view columnName="")
 Execute a user-defined reduce operation on the values of a column.
 
RResultPtr< RCutFlowReportReport ()
 Gather filtering statistics.
 
RResultPtr< RInterface< RLoopManager > > Snapshot (std::string_view treename, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options=RSnapshotOptions())
 Save selected columns to disk, in a new TTree treename in file filename.
 
RResultPtr< RInterface< RLoopManager > > Snapshot (std::string_view treename, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options=RSnapshotOptions())
 Save selected columns to disk, in a new TTree treename in file filename.
 
RResultPtr< RInterface< RLoopManager > > Snapshot (std::string_view treename, std::string_view filename, std::initializer_list< std::string > columnList, const RSnapshotOptions &options=RSnapshotOptions())
 Save selected columns to disk, in a new TTree treename in file filename.
 
RResultPtr< RInterface< RLoopManager > > Snapshot (std::string_view treename, std::string_view filename, std::string_view columnNameRegexp="", const RSnapshotOptions &options=RSnapshotOptions())
 Save selected columns to disk, in a new TTree treename in file filename.
 
RResultPtr< TStatisticStats (std::string_view value, std::string_view weight)
 Return a TStatistic object, filled once per event (lazy action)
 
RResultPtr< TStatisticStats (std::string_view value="")
 Return a TStatistic object, filled once per event (lazy action)
 
RResultPtr< doubleStdDev (std::string_view columnName="")
 Return the unbiased standard deviation of processed column values (lazy action)
 
RResultPtr< RDFDetail::SumReturnType_t< T > > Sum (std::string_view columnName="", const RDFDetail::SumReturnType_t< T > &initValue=RDFDetail::SumReturnType_t< T >{})
 Return the sum of processed column values (lazy action)
 
RResultPtr< COLL > Take (std::string_view column="")
 Return a collection of values of a column (lazy action, returns a std::vector by default)
 

Additional Inherited Members

- Protected Member Functions inherited from ROOT::RDF::RInterface< RDFDetail::RLoopManager >
 RInterface (const std::shared_ptr< RDFDetail::RLoopManager > &proxied, RLoopManager &lm, const RDFInternal::RBookedDefines &columns, RDataSource *ds)
 
void CheckAndFillDSColumns (ColumnNames_t validCols, TTraits::TypeList< ColumnTypes... > typeList)
 
RLoopManagerGetLoopManager () const
 
const std::shared_ptr< RDFDetail::RLoopManager > & GetProxiedPtr () const
 
ColumnNames_t GetValidatedColumnNames (const unsigned int nColumns, const ColumnNames_t &columns)
 Prepare the call to the GetValidatedColumnNames routine, making sure that GetBranchNames, which is expensive in terms of runtime, is called at most once.
 

#include <ROOT/RDataFrame.hxx>

Inheritance diagram for ROOT::RDataFrame:
[legend]

Member Typedef Documentation

◆ ColumnNames_t

Constructor & Destructor Documentation

◆ RDataFrame() [1/6]

ROOT::RDataFrame::RDataFrame ( std::string_view  treeName,
std::string_view  filenameglob,
const ColumnNames_t defaultBranches = {} 
)

Build the dataframe.

Parameters
[in]treeNameName of the tree contained in the directory
[in]filenameglobTDirectory where the tree is stored, e.g. a TFile.
[in]defaultBranchesCollection of default branches.

The filename globbing supports the same type of expressions as TChain::Add(). The default branches are looked at in case no branch is specified in the booking of actions or transformations. See RInterface for the documentation of the methods available.

Definition at line 1036 of file RDataFrame.cxx.

◆ RDataFrame() [2/6]

ROOT::RDataFrame::RDataFrame ( std::string_view  treeName,
const std::vector< std::string > &  fileglobs,
const ColumnNames_t defaultBranches = {} 
)

Build the dataframe.

Parameters
[in]treeNameName of the tree contained in the directory
[in]fileglobsCollection of file names of filename globs
[in]defaultBranchesCollection of default branches.

The filename globbing supports the same type of expressions as TChain::Add(). The default branches are looked at in case no branch is specified in the booking of actions or transformations. See RInterface for the documentation of the methods available.

Definition at line 1055 of file RDataFrame.cxx.

◆ RDataFrame() [3/6]

ROOT::RDataFrame::RDataFrame ( std::string_view  treeName,
::TDirectory dirPtr,
const ColumnNames_t defaultBranches = {} 
)

◆ RDataFrame() [4/6]

ROOT::RDataFrame::RDataFrame ( TTree tree,
const ColumnNames_t defaultBranches = {} 
)

Build the dataframe.

Parameters
[in]treeThe tree or chain to be studied.
[in]defaultBranchesCollection of default column names to fall back to when none is specified.

The default branches are looked at in case no branch is specified in the booking of actions or transformations. See RInterface for the documentation of the methods available.

Definition at line 1074 of file RDataFrame.cxx.

◆ RDataFrame() [5/6]

ROOT::RDataFrame::RDataFrame ( ULong64_t  numEntries)

Build a dataframe that generates numEntries entries.

Parameters
[in]numEntriesThe number of entries to generate.

An empty-source dataframe constructed with a number of entries will generate those entries on the fly when some action is triggered, and it will do so for all the previously-defined temporary branches. See RInterface for the documentation of the methods available.

Definition at line 1087 of file RDataFrame.cxx.

◆ RDataFrame() [6/6]

ROOT::RDataFrame::RDataFrame ( std::unique_ptr< ROOT::RDF::RDataSource ds,
const ColumnNames_t defaultBranches = {} 
)

Build dataframe associated to datasource.

Parameters
[in]dsThe data-source object.
[in]defaultBranchesCollection of default column names to fall back to when none is specified.

A dataframe associated to a datasource will query it to access column values. See RInterface for the documentation of the methods available.

Definition at line 1100 of file RDataFrame.cxx.

Libraries for ROOT::RDataFrame:

The documentation for this class was generated from the following files: