Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RDataFrame.cxx
Go to the documentation of this file.
1// Author: Enrico Guiraud, Danilo Piparo CERN 12/2016
2
3/*************************************************************************
4 * Copyright (C) 1995-2018, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11#include <algorithm>
12#include <stdexcept>
13
14#include "ROOT/RDataFrame.hxx"
15#include "ROOT/RDataSource.hxx"
16#include "TChain.h"
17#include "TDirectory.h"
18
19// clang-format off
20/**
21* \class ROOT::RDataFrame
22* \ingroup dataframe
23* \brief ROOT's RDataFrame offers a high level interface for analyses of data stored in `TTree`s, CSV's and other data formats.
24
25In addition, multi-threading and other low-level optimisations allow users to exploit all the resources available
26on their machines completely transparently.<br>
27Skip to the [class reference](#reference) or keep reading for the user guide.
28
29In a nutshell:
30~~~{.cpp}
31ROOT::EnableImplicitMT(); // Tell ROOT you want to go parallel
32ROOT::RDataFrame d("myTree", "file_*.root"); // Interface to TTree and TChain
33auto myHisto = d.Histo1D("Branch_A"); // This books the (lazy) filling of a histogram
34myHisto->Draw(); // Event loop is run here, upon first access to a result
35~~~
36
37Calculations are expressed in terms of a type-safe *functional chain of actions and transformations*, `RDataFrame` takes
38care of their execution. The implementation automatically puts in place several low level optimisations such as
39multi-thread parallelisation and caching.
40
41\htmlonly
42<a href="https://doi.org/10.5281/zenodo.260230"><img src="https://zenodo.org/badge/DOI/10.5281/zenodo.260230.svg"
43alt="DOI"></a>
44\endhtmlonly
45
46## For the impatient user
47You can directly see RDataFrame in action in our [tutorials](https://root.cern.ch/doc/master/group__tutorial__dataframe.html), in C++ or Python.
48
49## Table of Contents
50- [Cheat sheet](#cheatsheet)
51- [Introduction](#introduction)
52- [Crash course](#crash-course)
53- [Working with collections](#collections)
54- [Efficient analysis in Python](#python)
55- [Distributed execution in Python](#distrdf)
56- [Transformations](#transformations) -- manipulating data
57- [Actions](#actions) -- getting results
58- [Performance tips and parallel execution](#parallel-execution) -- how to use it and common pitfalls
59- [More features](#more-features)
60- [Class reference](#reference) -- most methods are implemented in the [RInterface](https://root.cern/doc/master/classROOT_1_1RDF_1_1RInterface.html) base class
61
62## <a name="cheatsheet"></a>Cheat sheet
63These are the operations which can be performed with RDataFrame.
64
65### Transformations
66Transformations are a way to manipulate the data.
67
68| **Transformation** | **Description** |
69|------------------|--------------------|
70| [Define](classROOT_1_1RDF_1_1RInterface.html#a7d48eb23b4378e99ebccb35e94ad025a) | 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 `pt`s of "good" muons). |
71| [DefineSlot](classROOT_1_1RDF_1_1RInterface.html#acaacf727b8a41d27c6bb4513348ac892) | 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`. |
72| [DefineSlotEntry](classROOT_1_1RDF_1_1RInterface.html#a4f17074d5771916e3df18f8458186de7) | 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. |
73| [Filter](classROOT_1_1RDF_1_1RInterface.html#a70284a3bedc72b19610aaa91b5007ebd) | Filter rows based on user-defined conditions. |
74| [Range](classROOT_1_1RDF_1_1RInterface.html#a1b36b7868831de2375e061bb06cfc225) | Filter rows based on entry number (single-thread only). |
75
76### Actions
77Actions aggregate data into a result. Each one is described in more detail in the reference guide.
78
79In 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.
80
81Lazy actions only trigger the event loop when one of the results is accessed for the first time, making it easy to
82produce many different results in one event loop. Instant actions trigger the event loop instantly.
83
84
85| **Lazy action** | **Description** |
86|------------------|-----------------|
87| [Aggregate](classROOT_1_1RDF_1_1RInterface.html#ae540b00addc441f9b504cbae0ef0a24d) | Execute a user-defined accumulation operation on the processed column values. |
88| [Book](https://root.cern/doc/master/classROOT_1_1RDF_1_1RInterface.html#a9ed8313806398d106bfc2390301c0408) | Book execution of a custom action using a user-defined helper object. |
89| [Cache](classROOT_1_1RDF_1_1RInterface.html#aaaa0a7bb8eb21315d8daa08c3e25f6c9) | 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). |
90| [Count](classROOT_1_1RDF_1_1RInterface.html#a37f9e00c2ece7f53fae50b740adc1456) | Return the number of events processed. Useful e.g. to get a quick count of the number of events passing a Filter. |
91| [Display](classROOT_1_1RDF_1_1RInterface.html#aee68f4411f16f00a1d46eccb6d296f01) | Provides a printable representation of the dataset contents. The method returns a [RDisplay](classROOT_1_1RDF_1_1RDisplay.html) instance which can be queried to get a compressed tabular representation on the standard output or a complete representation as a string. |
92| [Fill](classROOT_1_1RDF_1_1RInterface.html#a0cac4d08297c23d16de81ff25545440a) | Fill a user-defined object with the values of the specified columns, as if by calling `Obj.Fill(col1, col2, ...). |
93| [Graph](classROOT_1_1RDF_1_1RInterface.html#a804b466ebdbddef5c7e3400cc6b89301) | 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. |
94| [Histo{1D,2D,3D}](classROOT_1_1RDF_1_1RInterface.html#a247ca3aeb7ce5b95015b7fae72983055) | Fill a {one,two,three}-dimensional histogram with the processed column values. |
95| [Max](classROOT_1_1RDF_1_1RInterface.html#a057179b1e77599466a0b02200d5cd8c3) | Return the maximum of processed column values. If the type of the column is inferred, the return type is `double`, the type of the column otherwise.|
96| [Mean](classROOT_1_1RDF_1_1RInterface.html#ade6b020284f2f4fe9d3b09246b5f376a) | Return the mean of processed column values.|
97| [Min](classROOT_1_1RDF_1_1RInterface.html#a7005702189e601972b6d19ecebcdc80c) | Return the minimum of processed column values. If the type of the column is inferred, the return type is `double`, the type of the column otherwise.|
98| [Profile{1D,2D}](classROOT_1_1RDF_1_1RInterface.html#a8ef7dc16b0e9f7bc9cfbe2d9e5de0cef) | Fill a {one,two}-dimensional profile with the column values that passed all filters. |
99| [Reduce](classROOT_1_1RDF_1_1RInterface.html#a118e723ae29834df8f2a992ded347354) | Reduce (e.g. sum, merge) entries using the function (lambda, functor...) passed as argument. The function must have signature `T(T,T)` where `T` is 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. |
100| [Report](classROOT_1_1RDF_1_1RInterface.html#a94f322531dcb25beb8f53a602e5d6332) | 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. |
101| [Stats](https://root.cern/doc/master/classROOT_1_1RDF_1_1RInterface.html#a9e8fafb75abfa4faed4da18dcde01568) | Return a TStatistic object filled with the input columns. |
102| [StdDev](classROOT_1_1RDF_1_1RInterface.html#a482c4e4f81fe1e421c016f89cd281572) | Return the unbiased standard deviation of the processed column values. |
103| [Sum](classROOT_1_1RDF_1_1RInterface.html#a61d03407459120df6749af43ed506891) | Return the sum of the values in the column. If the type of the column is inferred, the return type is `double`, the type of the column otherwise. |
104| [Take](classROOT_1_1RDF_1_1RInterface.html#a4fd694773a2931b6b07737ddcd1e73b4) | 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 a `ROOT::VecOps::RVec<T>` to guarantee the lifetime of the data involved. |
105
106| **Instant action** | **Description** |
107|---------------------|-----------------|
108| [Foreach](classROOT_1_1RDF_1_1RInterface.html#ad2822a7ccb8a9afdf3e5b2ea321886ca) | 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. |
109| [ForeachSlot](classROOT_1_1RDF_1_1RInterface.html#a3650ca30aae1ccd0d92bf3d680314129) | 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`. |
110| [Snapshot](classROOT_1_1RDF_1_1RInterface.html#a233b7723e498967f4340705d2c4db7f8) | 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.|
111
112
113### Other Operations
114
115| **Operation** | **Description** |
116|---------------------|-----------------|
117| [Alias](classROOT_1_1RDF_1_1RInterface.html#a31ca327e4a192dcc05a4aac240e1a725) | Introduce an alias for a particular column name. |
118| [GetColumnNames](classROOT_1_1RDF_1_1RInterface.html#a951fe60b74d3a9fda37df59fd1dac186) | Get the names of all the available columns of the dataset. |
119| [GetDefinedColumnNames](classROOT_1_1RDF_1_1RInterface.html#ad5c3fab8155aae8f614735df68430c58) | Get the names of all the defined columns |
120| [GetColumnType](classROOT_1_1RDF_1_1RInterface.html#ad3ccd813d9fed014ae6a080411c5b5a8) | Return the type of a given column as a string. |
121| [GetColumnTypeNamesList](classROOT_1_1RDF_1_1RInterface.html#a951fe60b74d3a9fda37df59fd1dac186) | Return the list of type names of columns in the dataset. |
122| [GetFilterNames](classROOT_1_1RDF_1_1RInterface.html#a25026681111897058299161a70ad9bb2) | 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. |
123| [SaveGraph](https://root.cern/doc/master/namespaceROOT_1_1RDF.html#a366b19a07428c69801020e2edba117dd) | Store the computation graph of an RDataFrame in graphviz format for easy inspection. |
124| [GetNRuns](classROOT_1_1RDF_1_1RInterface.html#adfb0562a9f7732c3afb123aefa07e0df) | Get the number of event loops run by this RDataFrame instance. |
125
126
127## <a name="introduction"></a>Introduction
128Users define their analysis as a sequence of operations to be performed on the data-frame object; the framework
129takes care of the management of the loop over entries as well as low-level details such as I/O and parallelisation.
130`RDataFrame` provides methods to perform most common operations required by ROOT analyses;
131at the same time, users can just as easily specify custom code that will be executed in the event loop.
132
133`RDataFrame` is built with a *modular* and *flexible* workflow in mind, summarised as follows:
134
1351. **build a data-frame** object by specifying your data-set
1362. **apply a series of transformations** to your data
137 1. **filter** (e.g. apply some cuts) or
138 2. **define** a new column (e.g. the result of an expensive computation on columns)
1393. **apply actions** to the transformed data to produce results (e.g. fill a histogram)
140
141Make 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.
142
143The following table shows how analyses based on `TTreeReader` and `TTree::Draw` translate to `RDataFrame`. Follow the
144[crash course](#crash-course) to discover more idiomatic and flexible ways to express analyses with `RDataFrame`.
145<table>
146<tr>
147 <td>
148 <b>TTreeReader</b>
149 </td>
150 <td>
151 <b>ROOT::RDataFrame</b>
152 </td>
153</tr>
154<tr>
155 <td>
156~~~{.cpp}
157TTreeReader reader("myTree", file);
158TTreeReaderValue<A_t> a(reader, "A");
159TTreeReaderValue<B_t> b(reader, "B");
160TTreeReaderValue<C_t> c(reader, "C");
161while(reader.Next()) {
162 if(IsGoodEvent(*a, *b, *c))
163 DoStuff(*a, *b, *c);
164}
165~~~
166 </td>
167 <td>
168~~~{.cpp}
169ROOT::RDataFrame d("myTree", file, {"A", "B", "C"});
170d.Filter(IsGoodEvent).Foreach(DoStuff);
171~~~
172 </td>
173</tr>
174<tr>
175 <td>
176 <b>TTree::Draw</b>
177 </td>
178 <td>
179 <b>ROOT::RDataFrame</b>
180 </td>
181</tr>
182<tr>
183 <td>
184~~~{.cpp}
185auto t = file->Get<TTree>("myTree");
186t->Draw("x", "y > 2");
187~~~
188 </td>
189 <td>
190~~~{.cpp}
191ROOT::RDataFrame d("myTree", file);
192auto h = d.Filter("y > 2").Histo1D("x");
193~~~
194 </td>
195</tr>
196</table>
197
198## <a name="crash-course"></a> Crash course
199All snippets of code presented in the crash course can be executed in the ROOT interpreter. Simply precede them with
200~~~{.cpp}
201using namespace ROOT; // RDataFrame's namespace
202~~~
203which is omitted for brevity. The terms "column" and "branch" are used interchangeably.
204
205### Creating a RDataFrame
206RDataFrame's constructor is where the user specifies the dataset and, optionally, a default set of columns that
207operations should work with. Here are the most common methods to construct a RDataFrame object:
208~~~{.cpp}
209// single file -- all ctors are equivalent
210TFile *f = TFile::Open("file.root");
211TTree *t = f.Get<TTree>("treeName");
212
213RDataFrame d1("treeName", "file.root");
214RDataFrame d2("treeName", f); // same as TTreeReader
215RDataFrame d3(*t);
216
217// multiple files -- all ctors are equivalent
218TChain chain("myTree");
219chain.Add("file1.root");
220chain.Add("file2.root");
221
222RDataFrame d4("myTree", {"file1.root", "file2.root"});
223std::vector<std::string> files = {"file1.root", "file2.root"};
224RDataFrame d5("myTree", files);
225RDataFrame d6("myTree", "file*.root"); // the glob is passed as-is to TChain's constructor
226RDataFrame d7(chain);
227~~~
228Additionally, users can construct a RDataFrame with no data source by passing an integer number. This is the number of rows that
229will be generated by this RDataFrame.
230~~~{.cpp}
231RDataFrame d(10); // a RDF with 10 entries (and no columns/branches, for now)
232d.Foreach([] { static int i = 0; std::cout << i++ << std::endl; }); // silly example usage: count to ten
233~~~
234This 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).
235
236For data sources other than TTrees and TChains, RDataFrame objects are constructed using ad-hoc factory functions (see e.g. MakeCsvDataFrame(), MakeSqliteDataFrame(), MakeArrowDataFrame()):
237
238~~~{.cpp}
239auto df = ROOT::RDF::MakeCsvDataFrame("input.csv");
240// use df as usual
241~~~
242
243### Filling a histogram
244Let's now tackle a very common task, filling a histogram:
245~~~{.cpp}
246// Fill a TH1D with the "MET" branch
247RDataFrame d("myTree", "file.root");
248auto h = d.Histo1D("MET");
249h->Draw();
250~~~
251The first line creates a `RDataFrame` associated to the `TTree` "myTree". This tree has a branch named "MET".
252
253Histo1D() is an *action*; it returns a smart pointer (a RResultPtr, to be precise) to a TH1D histogram filled
254with the `MET` of all events. If the quantity stored in the branch is a collection (e.g. a vector or an array), the
255histogram is filled with all vector elements for each event.
256
257You can use the objects returned by actions as if they were pointers to the desired results. There are many other
258possible [actions](#cheatsheet), and all their results are wrapped in smart pointers; we'll see why in a minute.
259
260### Applying a filter
261Let'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:
262~~~{.cpp}
263RDataFrame d("myTree", "file.root");
264auto c = d.Filter("MET > 4.").Count(); // computations booked, not run
265std::cout << *c << std::endl; // computations run here, upon first access to the result
266~~~
267The filter string (which must contain a valid C++ expression) is applied to the specified branches for each event;
268the name and types of the columns are inferred automatically. The string expression is required to return a `bool`
269which signals whether the event passes the filter (`true`) or not (`false`).
270
271You can think of your data as "flowing" through the chain of calls, being transformed, filtered and finally used to
272perform actions. Multiple `Filter` calls can be chained one after another.
273
274Using string filters is nice for simple things, but they are limited to specifying the equivalent of a single return
275statement or the body of a lambda, so it's cumbersome to use strings with more complex filters. They also add a small
276runtime overhead, as ROOT needs to just-in-time compile the string into C++ code. When more freedom is required or
277runtime performance is very important, a C++ callable can be specified instead (a lambda in the following snippet,
278but it can be any kind of function or even a functor class), together with a list of branch names.
279This snippet is analogous to the one above:
280~~~{.cpp}
281RDataFrame d("myTree", "file.root");
282auto metCut = [](double x) { return x > 4.; }; // a C++11 lambda function checking "x > 4"
283auto c = d.Filter(metCut, {"MET"}).Count();
284std::cout << *c << std::endl;
285~~~
286
287An example of a more complex filter expressed as a string containing C++ code is shown below
288
289~~~{.cpp}
290RDataFrame d("myTree", "file.root");
291auto df = d.Define("p", "std::array<double, 4> p{px, py, pz}; return p;")
292 .Filter("double p2 = 0.0; for (auto&& x : p) p2 += x*x; return sqrt(p2) < 10.0;");
293~~~
294
295The code snippet above defines a column `p` that is a fixed-size array using the component column names and then
296filters on its magnitude by looping over its elements. It must be noted that the usage of strings to define columns
297like the one above is a major advantage when using PyROOT. However, only constants and data coming from other columns
298in the dataset can be involved in the code passed as a string. Local variables and functions cannot be used, since
299the interpreter will not know how to find them. When capturing local state is necessary, a C++ callable can be used.
300
301More information on filters and how to use them to automatically generate cutflow reports can be found [below](#Filters).
302
303### Defining custom columns
304Let's now consider the case in which "myTree" contains two quantities "x" and "y", but our analysis relies on a derived
305quantity `z = sqrt(x*x + y*y)`. Using the `Define` transformation, we can create a new column in the data-set containing
306the variable "z":
307~~~{.cpp}
308RDataFrame d("myTree", "file.root");
309auto sqrtSum = [](double x, double y) { return sqrt(x*x + y*y); };
310auto zMean = d.Define("z", sqrtSum, {"x","y"}).Mean("z");
311std::cout << *zMean << std::endl;
312~~~
313`Define` creates the variable "z" by applying `sqrtSum` to "x" and "y". Later in the chain of calls we refer to
314variables created with `Define` as if they were actual tree branches/columns, but they are evaluated on demand, at most
315once per event. As with filters, `Define` calls can be chained with other transformations to create multiple custom
316columns. `Define` and `Filter` transformations can be concatenated and intermixed at will.
317
318As with filters, it is possible to specify new columns as string expressions. This snippet is analogous to the one above:
319~~~{.cpp}
320RDataFrame d("myTree", "file.root");
321auto zMean = d.Define("z", "sqrt(x*x + y*y)").Mean("z");
322std::cout << *zMean << std::endl;
323~~~
324Again the names of the branches used in the expression and their types are inferred automatically. The string must be
325valid C++ and is just-in-time compiled by the ROOT interpreter, cling -- the process has a small runtime overhead.
326
327Previously, when showing the different ways a RDataFrame can be created, we showed a constructor that only takes a
328number of entries a parameter. In the following example we show how to combine such an "empty" `RDataFrame` with `Define`
329transformations to create a data-set on the fly. We then save the generated data on disk using the `Snapshot` action.
330~~~{.cpp}
331RDataFrame d(100); // a RDF that will generate 100 entries (currently empty)
332int x = -1;
333auto d_with_columns = d.Define("x", [&x] { return ++x; })
334 .Define("xx", [&x] { return x*x; });
335d_with_columns.Snapshot("myNewTree", "newfile.root");
336~~~
337This example is slightly more advanced than what we have seen so far: for starters, it makes use of lambda captures (a
338simple way to make external variables available inside the body of C++ lambdas) to act on the same variable `x` from
339both `Define` transformations. Secondly we have *stored* the transformed data-frame in a variable. This is always
340possible: at each point of the transformation chain, users can store the status of the data-frame for further use (more
341on this [below](#callgraphs)).
342
343You can read more about defining new columns [here](#custom-columns).
344
345\image html RDF_Graph.png "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."
346
347### Running on a range of entries
348It is sometimes necessary to limit the processing of the dataset to a range of entries. For this reason, the RDataFrame
349offers the concept of ranges as a node of the RDataFrame chain of transformations; this means that filters, columns and
350actions can be concatenated to and intermixed with `Range`s. If a range is specified after a filter, the range will act
351exclusively on the entries passing the filter -- it will not even count the other entries! The same goes for a `Range`
352hanging from another `Range`. Here are some commented examples:
353~~~{.cpp}
354RDataFrame d("myTree", "file.root");
355// Here we store a data-frame that loops over only the first 30 entries in a variable
356auto d30 = d.Range(30);
357// This is how you pick all entries from 15 onwards
358auto d15on = d.Range(15, 0);
359// We can specify a stride too, in this case we pick an event every 3
360auto d15each3 = d.Range(0, 15, 3);
361~~~
362Note that ranges are not available when multi-threading is enabled. More information on ranges is available
363[here](#ranges).
364
365### Executing multiple actions in the same event loop
366As a final example let us apply two different cuts on branch "MET" and fill two different histograms with the "pt\_v" of
367the filtered events.
368By now, you should be able to easily understand what's happening:
369~~~{.cpp}
370RDataFrame d("treeName", "file.root");
371auto h1 = d.Filter("MET > 10").Histo1D("pt_v");
372auto h2 = d.Histo1D("pt_v");
373h1->Draw(); // event loop is run once here
374h2->Draw("SAME"); // no need to run the event loop again
375~~~
376`RDataFrame` executes all above actions by **running the event-loop only once**. The trick is that actions are not
377executed at the moment they are called, but they are **lazy**, i.e. delayed until the moment one of their results is
378accessed through the smart pointer. At that time, the event loop is triggered and *all* results are produced
379simultaneously.
380
381It is therefore good practice to declare all your transformations and actions *before* accessing their results, allowing
382`RDataFrame` to run the loop once and produce all results in one go.
383
384### Going parallel
385Let's say we would like to run the previous examples in parallel on several cores, dividing events fairly between cores.
386The only modification required to the snippets would be the addition of this line *before* constructing the main
387data-frame object:
388~~~{.cpp}
389ROOT::EnableImplicitMT();
390~~~
391Simple as that. More details are given [below](#parallel-execution).
392
393## <a name="collections"></a> Working with collections
394
395RDataFrame reads collections as the special type RVec (e.g. a branch containing an array of floating point numbers can
396be read as a `RVec<float>`). C-style arrays (with variable or static size), `std::vector`s and most other collection
397types can be read this way. When reading ROOT data, column values of type `RVec<T>` perform no copy of the underlying array.
398
399RVec is a container similar to `std::vector` but it offers a rich interface to operate on the array elements in a
400vectorised fashion, similar to Python's NumPy arrays.
401
402For example, to fill a histogram with the `pt` of selected particles for each event, Define() can be used to create
403a column that contains the desired array elements as follows:
404
405~~~{.cpp}
406# h is filled with all the elements of `good_pts`, for each event
407h = df.Define("good_pts", "pt[pt > 0]").Histo1D("good_pts")
408~~~
409
410Learn more on [RVec](https://root.cern/doc/master/classROOT_1_1VecOps_1_1RVec.html).
411
412## <a name="python"></a>Efficient analysis in Python
413
414You can use `RDataFrame` in Python due to the dynamic C++/Python translation of PyROOT. In general, the interface
415is the same as for C++, a simple example follows.
416
417~~~{.python}
418df = ROOT.RDataFrame("myTree", "myFile.root")
419sum = df.Filter("x > 10").Sum("y")
420print(sum.GetValue())
421~~~
422
423### Simple usage of efficient C++ code in Python
424
425To perform more complex operations in the `RDataFrame` graph, e.g., in `Filter` and `Define` nodes, which don't
426fit into a simple expression string, you can just-in-time compile such functions directly in the Python script
427via the C++ interpreter cling. This approach has the advantage that you get the efficiency of compiled C++ code
428combined with the convenient workflow of a Python script. See the following snippet for an example of how to
429use a just-in-time-compiled C++ function from Python.
430
431~~~{.python}
432ROOT.gInterpreter.Declare("""
433bool myFilter(float x) {
434 return x > 10;
435}
436""")
437
438df = ROOT.RDataFrame("myTree", "myFile.root")
439sum = df.Filter("myFilter(x)").Sum("y")
440print(sum.GetValue())
441~~~
442
443To increase the performance even further, you can also pre-compile a C++ library with full code optimizations
444and load the function into the `RDataFrame` computation as follows.
445
446~~~{.python}
447ROOT.gSystem.Load("path/to/myLibrary.so") # Library with the myFilter function
448ROOT.gInterpreter.Declare('#include "myLibrary.h"') # Header with the definition of the myFilter function
449df = ROOT.RDataFrame("myTree", "myFile.root")
450sum = df.Filter("myFilter(x)").Sum("y")
451print(sum.GetValue())
452~~~
453
454### Just-in-time compilation of Python callables with numba
455
456ROOT also offers the option to compile Python callables with fundamental types and arrays thereof using numba and then
457using the function in `RDataFrame` from C++. The workflow requires the Python packages `numba` and `cffi`
458to be installed. See the following snippet for a simple example or the full tutorial [here](pyroot004__NumbaDeclare_8py.html).
459
460~~~{.python}
461@ROOT.Numba.Declare(["float"], "bool")
462def myFilter(x):
463 return x > 10
464
465df = ROOT.RDataFrame("myTree", "myFile.root")
466sum = df.Filter("Numba::myFilter(x)").Sum("y")
467print(sum.GetValue())
468~~~
469
470### Conversion to numpy arrays
471
472Eventually, you probably would like to inspect the content of the `RDataFrame` or process the data further
473with functionality from Python libraries. For this purpose, we provide the `AsNumpy` function, which is able
474to provide you the columns of your `RDataFrame` as numpy arrays in Python. See a brief introduction below or
475a full tutorial [here](df026__AsNumpyArrays_8py.html).
476
477~~~{.python}
478df = ROOT.RDataFrame("myTree", "myFile.root")
479cols = df.Filter("x > 10").AsNumpy(["x", "y"])
480print(cols["x"], cols["y"])
481~~~
482
483## <a name="distrdf"></a>Distributed execution in Python
484
485RDataFrame applications can be executed in parallel through distributed computing frameworks on a set of remote machines
486thanks to the Python package `ROOT.RDF.Experimental.Distributed`. This experimental, **Python-only** package allows to scale the
487optimized performance RDataFrame can achieve on a single machine to multiple nodes at the same time. It is designed so
488that different backends can be easily plugged in, currently supporting [Apache Spark](http://spark.apache.org/) and soon
489also [Dask](https://dask.org/). To make use of distributed RDataFrame, you only need to switch `ROOT.RDataFrame` with
490the backend-specific `RDataFrame` of your choice, for example:
491
492~~~{.py}
493import ROOT
494
495# Point RDataFrame calls to the Spark specific RDataFrame
496RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
497
498# It still accepts the same constructor arguments as traditional RDataFrame
499df = RDataFrame("mytree", "myfile.root")
500
501# Continue the application with the traditional RDataFrame API
502sum = df.Filter("x > 10").Sum("y")
503h = df.Histo1D("x")
504
505print(sum.GetValue())
506h.Draw()
507~~~
508
509The main goal of this package is to support running any RDataFrame application distributedly. Nonetheless, not all
510RDataFrame operations currently work with this package. The subset that is currently available is:
511- AsNumpy
512- Count
513- Define
514- Fill
515- Filter
516- Graph
517- Histo[1,2,3]D
518- Max
519- Mean
520- Min
521- Profile[1,2,3]D
522- Snapshot
523- Sum
524
525with support for more operations coming in the future. Data sources other than TTree and TChain (e.g. CSV, RNTuple) are
526currently not supported.
527
528### Connecting to a Spark cluster
529
530In order to distribute the RDataFrame workload, you can connect to a Spark cluster you have access to through the
531official [Spark API](https://spark.apache.org/docs/latest/rdd-programming-guide.html#initializing-spark), then hook the
532connection instance to the distributed `RDataFrame` object like so:
533
534~~~{.py}
535import pyspark
536import ROOT
537
538# Create a SparkContext object with the right configuration for your Spark cluster
539conf = SparkConf().setAppName(appName).setMaster(master)
540sc = SparkContext(conf=conf)
541
542# Point RDataFrame calls to the Spark specific RDataFrame
543RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
544
545# The Spark RDataFrame constructor accepts an optional "sparkcontext" parameter
546# and it will distribute the application to the connected cluster
547df = RDataFrame("mytree", "myfile.root", sparkcontext = sc)
548~~~
549
550If an instance of [SparkContext](https://spark.apache.org/docs/latest/api/python/reference/api/pyspark.SparkContext.html)
551is not provided, the default behaviour is to create one in the background for you.
552
553## <a name="transformations"></a>Transformations
554### <a name="Filters"></a> Filters
555A filter is created through a call to `Filter(f, columnList)` or `Filter(filterString)`. In the first overload, `f` can
556be a function, a lambda expression, a functor class, or any other callable object. It must return a `bool` signalling
557whether the event has passed the selection (`true`) or not (`false`). It should perform "read-only" operations on the
558columns, and should not have side-effects (e.g. modification of an external or static variable) to ensure correctness
559when implicit multi-threading is active. The second overload takes a string with a valid C++ expression in which column
560names are used as variable names (e.g. `Filter("x[0] + x[1] > 0")`). This is a convenience feature that comes with a
561certain runtime overhead: C++ code has to be generated on the fly from this expression before using it in the event
562loop. See the paragraph about "Just-in-time compilation" below for more information.
563
564`RDataFrame` only evaluates filters when necessary: if multiple filters are chained one after another, they are executed
565in order and the first one returning `false` causes the event to be discarded and triggers the processing of the next
566entry. If multiple actions or transformations depend on the same filter, that filter is not executed multiple times for
567each entry: after the first access it simply serves a cached result.
568
569#### <a name="named-filters-and-cutflow-reports"></a>Named filters and cutflow reports
570An optional string parameter `name` can be passed to the `Filter` method to create a **named filter**. Named filters
571work as usual, but also keep track of how many entries they accept and reject.
572
573Statistics are retrieved through a call to the `Report` method:
574
575- when `Report` is called on the main `RDataFrame` object, it returns a RResultPtr<RCutFlowReport> relative to all
576named filters declared up to that point
577- when called on a specific node (e.g. the result of a `Define` or `Filter`), it returns a RResultPtr<RCutFlowReport>
578relative all named filters in the section of the chain between the main `RDataFrame` and that node (included).
579
580Stats are stored in the same order as named filters have been added to the graph, and *refer to the latest event-loop*
581that has been run using the relevant `RDataFrame`.
582
583### <a name="ranges"></a>Ranges
584When `RDataFrame` is not being used in a multi-thread environment (i.e. no call to `EnableImplicitMT` was made),
585`Range` transformations are available. These act very much like filters but instead of basing their decision on
586a filter expression, they rely on `begin`,`end` and `stride` parameters.
587
588- `begin`: initial entry number considered for this range.
589- `end`: final entry number (excluded) considered for this range. 0 means that the range goes until the end of the dataset.
590- `stride`: process one entry of the [begin, end) range every `stride` entries. Must be strictly greater than 0.
591
592The actual number of entries processed downstream of a `Range` node will be `(end - begin)/stride` (or less if less
593entries than that are available).
594
595Note that ranges act "locally", not based on the global entry count: `Range(10,50)` means "skip the first 10 entries
596*that reach this node*, let the next 40 entries pass, then stop processing". If a range node hangs from a filter node,
597and the range has a `begin` parameter of 10, that means the range will skip the first 10 entries *that pass the
598preceding filter*.
599
600Ranges allow "early quitting": if all branches of execution of a functional graph reached their `end` value of
601processed entries, the event-loop is immediately interrupted. This is useful for debugging and quick data explorations.
602
603### <a name="custom-columns"></a> Custom columns
604Custom columns are created by invoking `Define(name, f, columnList)`. As usual, `f` can be any callable object
605(function, lambda expression, functor class...); it takes the values of the columns listed in `columnList` (a list of
606strings) as parameters, in the same order as they are listed in `columnList`. `f` must return the value that will be
607assigned to the temporary column.
608
609A new variable is created called `name`, accessible as if it was contained in the dataset from subsequent
610transformations/actions.
611
612Use cases include:
613- caching the results of complex calculations for easy and efficient multiple access
614- extraction of quantities of interest from complex objects
615- branch aliasing, i.e. changing the name of a branch
616
617An exception is thrown if the `name` of the new column/branch is already in use for another branch in the `TTree`.
618
619It is also possible to specify the quantity to be stored in the new temporary column as a C++ expression with the method
620`Define(name, expression)`. For example this invocation
621
622~~~{.cpp}
623df.Define("pt", "sqrt(px*px + py*py)");
624~~~
625
626will create a new column called "pt" the value of which is calculated starting from the columns px and py. The system
627builds a just-in-time compiled function starting from the expression after having deduced the list of necessary branches
628from the names of the variables specified by the user.
629
630#### Custom columns as function of slot and entry number
631
632It is possible to create custom columns also as a function of the processing slot and entry numbers. The methods that can
633be invoked are:
634- `DefineSlot(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, T1, T2, ...)`: the
635first parameter is the slot number which ranges from 0 to ROOT::GetThreadPoolSize() - 1.
636- `DefineSlotEntry(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, ULong64_t,
637T1, T2, ...)`: the first parameter is the slot number while the second one the number of the entry being processed.
638
639## <a name="actions"></a>Actions
640### Instant and lazy actions
641Actions can be **instant** or **lazy**. Instant actions are executed as soon as they are called, while lazy actions are
642executed whenever the object they return is accessed for the first time. As a rule of thumb, actions with a return value
643are lazy, the others are instant.
644
645## <a name="parallel-execution"></a>Performance tips and parallel execution
646As pointed out before in this document, RDataFrame can transparently perform multi-threaded event loops to speed up
647the execution of its actions. Users have to call ROOT::EnableImplicitMT() *before* constructing the `RDataFrame`
648object to indicate that it should take advantage of a pool of worker threads. **Each worker thread processes a distinct
649subset of entries**, and their partial results are merged before returning the final values to the user.
650There are no guarantees on the order in which threads will process the batches of entries.
651In particular, note that this means that, for multi-thread event loops, there is no
652guarantee on the order in which `Snapshot` will _write_ entries: they could be scrambled with respect to the input dataset.
653
654\warning By default, RDataFrame will use as many threads as the hardware supports, using up **all** the resources on
655a machine. This might be undesirable on shared computing resources such as a batch cluster. Therefore, when running on shared computing resources, use
656```
657ROOT::EnableImplicitMT(i)
658```
659replacing `i` with the number of CPUs/slots that were allocated for this job.
660
661### Thread-safety of user-defined expressions
662RDataFrame operations such as `Histo1D` or `Snapshot` are guaranteed to work correctly in multi-thread event loops.
663User-defined expressions, such as strings or lambdas passed to `Filter`, `Define`, `Foreach`, `Reduce` or `Aggregate`
664will have to be thread-safe, i.e. it should be possible to call them concurrently from different threads.
665
666Note that simple `Filter` and `Define` transformations will inherently satisfy this requirement: `Filter`/`Define`
667expressions will often be *pure* in the functional programming sense (no side-effects, no dependency on external state),
668which eliminates all risks of race conditions.
669
670In order to facilitate writing of thread-safe operations, some RDataFrame features such as `Foreach`, `Define` or `OnPartialResult`
671offer thread-aware counterparts (`ForeachSlot`, `DefineSlot`, `OnPartialResultSlot`): their only difference is that they
672will pass an extra `slot` argument (an unsigned integer) to the user-defined expression. When calling user-defined code
673concurrently, `RDataFrame` guarantees that different threads will employ different values of the `slot` parameter,
674where `slot` will be a number between 0 and `ROOT::GetThreadPoolSize() - 1`.
675In other words, within a slot, computation runs sequentially and events are processed sequentially.
676Note that the same slot might be associated to different threads over the course of a single event loop, but two threads
677will never receive the same slot at the same time.
678This extra parameter might facilitate writing safe parallel code by having each thread write/modify a different
679*processing slot*, e.g. a different element of a list. See [here](#generic-actions) for an example usage of `ForeachSlot`.
680
681### Parallel execution of multiple `RDataFrame` event loops
682A complex analysis may require multiple separate `RDatFrame` computation graphs to produce all desired results. This poses the challenge that the
683event loops of each computation graph can be parallelized, but the different loops run sequentially, one after the other.
684On many-core architectures it might be desirable to run different event loops concurrently to improve resource usage.
685ROOT::RDF::RunGraphs() allows running multiple `RDataFrame` event loops concurrently:
686~~~{.cpp}
687ROOT::EnableImplicitMT();
688ROOT::RDataFrame df1("tree1", "f1.root");
689ROOT::RDataFrame df2("tree2", "f2.root");
690auto histo1 = df1.Histo1D("x");
691auto histo2 = df2.Histo1D("y");
692
693// just accessing result pointers, the event loops of separate RDataFrames run one after the other
694histo1->Draw(); // runs first multi-thread event loop
695histo2->Draw(); // runs second multi-thread event loop
696
697// with ROOT::RDF::RunGraphs, event loops for separate computation graphs can run concurrently
698ROOT::RDF::RunGraphs({histo1, histo2});
699~~~
700
701### Performance profiling of RDataFrame applications
702
703To obtain the maximum performance out of RDataFrame, make sure to avoid just-in-time compiled versions of transformations and actions if at all possible.
704For instance, `Filter("x > 0")` requires just-in-time compilation of the corresponding C++ logic, while the equivalent `Filter([] { return x > 0; }, {"x"})` does not.
705Similarly, `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
706should be preferred for performance-critical applications.
707
708Python applications cannot easily specify template parameters or pass C++ callables to RDataFrame.
709See [Efficient analysis in Python](#python) for possible ways to speed up hot paths in this case.
710
711Just-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*
712before 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.
713
714Also 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:
715~~~{.cpp}
716auto verbosity = ROOT::Experimental::RLogScopedVerbosity(ROOT::Detail::RDF::RDFLogChannel(), ROOT::Experimental::ELogLevel.kInfo);
717~~~
718
719### Memory usage
720
721There 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.
722Secondly, 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.
723
724## <a name="more-features"></a>More features
725Here is a list of the most important features that have been omitted in the "Crash course" for brevity.
726You don't need to read all these to start using `RDataFrame`, but they are useful to save typing time and runtime.
727
728### Programmatically get the list of column names
729The `GetColumnsNames()` method returns the list of valid column names for the dataset:
730~~~{.cpp}
731RDataFrame d("myTree", "file.root");
732std::vector<std::string> colNames = d.GetColumnNames();
733~~~
734
735### Callbacks
736It's possible to schedule execution of arbitrary functions (callbacks) during the event loop.
737Callbacks can be used e.g. to inspect partial results of the analysis while the event loop is running,
738drawing a partially-filled histogram every time a certain number of new entries is processed, or event
739displaying a progress bar while the event loop runs.
740
741For example one can draw an up-to-date version of a result histogram every 100 entries like this:
742~~~{.cpp}
743auto h = tdf.Histo1D("x");
744TCanvas c("c","x hist");
745h.OnPartialResult(100, [&c](TH1D &h_) { c.cd(); h_.Draw(); c.Update(); });
746h->Draw(); // event loop runs here, this `Draw` is executed after the event loop is finished
747~~~
748
749Callbacks are registered to a RResultPtr and must be callables that takes a reference to the result type as argument
750and return nothing. RDataFrame will invoke registered callbacks passing partial action results as arguments to them
751(e.g. a histogram filled with a part of the selected events).
752
753Read more on RResultPtr::OnPartialResult().
754
755### Default branch lists
756When constructing a `RDataFrame` object, it is possible to specify a **default column list** for your analysis, in the
757usual form of a list of strings representing branch/column names. The default column list will be used as a fallback
758whenever a list specific to the transformation/action is not present. RDataFrame will take as many of these columns as
759needed, ignoring trailing extra names if present.
760~~~{.cpp}
761// use "b1" and "b2" as default branches
762RDataFrame d1("myTree", "file.root", {"b1","b2"});
763auto h = d1.Filter([](int b1, int b2) { return b1 > b2; }) // will act on "b1" and "b2"
764 .Histo1D(); // will act on "b1"
765
766// just one default branch this time
767RDataFrame d2("myTree", "file.root", {"b1"});
768auto min = d2.Filter([](double b2) { return b2 > 0; }, {"b2"}) // we can still specify non-default branch lists
769 .Min(); // returns the minimum value of "b1" for the filtered entries
770~~~
771
772### <a name="ImplicitColumns"></a> Implicit Columns
773Every instance of `RDataFrame` is created with two special columns called `rdfentry_` and `rdfslot_`. The `rdfentry_`
774column is an unsigned 64-bit integer holding the current entry number while `rdfslot_` is an unsigned 32-bit integer
775holding the index of the current data processing slot.
776For backwards compatibility reasons, the names `tdfentry_` and `tdfslot_` are also accepted.
777These columns are not considered by operations such as [Cache](classROOT_1_1RDF_1_1RInterface.html#aaaa0a7bb8eb21315d8daa08c3e25f6c9)
778or [Snapshot](classROOT_1_1RDF_1_1RInterface.html#a233b7723e498967f4340705d2c4db7f8). The _cached_ or _snapshot_ data frame
779provides "its own" values for these columns which do not necessarily correspond to the ones of the mother data frame. This is
780most notably the case where filters are used before deriving a cached/persistified dataframe.
781
782Note that in multi-thread event loops the values of `rdfentry_` _do not_ correspond to what would be the entry numbers
783of a TChain constructed over the same set of ROOT files, as the entries are processed in an unspecified order.
784
785### Just-in-time compilation: branch type inference and explicit declaration of branch types
786C++ is a statically typed language: all types must be known at compile-time. This includes the types of the `TTree`
787branches we want to work on. For filters, temporary columns and some of the actions, **branch types are deduced from the
788signature** of the relevant filter function/temporary column expression/action function:
789~~~{.cpp}
790// here b1 is deduced to be `int` and b2 to be `double`
791dataFrame.Filter([](int x, double y) { return x > 0 && y < 0.; }, {"b1", "b2"});
792~~~
793If we specify an incorrect type for one of the branches, an exception with an informative message will be thrown at
794runtime, when the branch value is actually read from the `TTree`: `RDataFrame` detects type mismatches. The same would
795happen if we swapped the order of "b1" and "b2" in the branch list passed to `Filter`.
796
797Certain actions, on the other hand, do not take a function as argument (e.g. `Histo1D`), so we cannot deduce the type of
798the branch at compile-time. In this case **`RDataFrame` infers the type of the branch** from the `TTree` itself. This
799is why we never needed to specify the branch types for all actions in the above snippets.
800
801When the branch type is not a common one such as `int`, `double`, `char` or `float` it is nonetheless good practice to
802specify it as a template parameter to the action itself, like this:
803~~~{.cpp}
804dataFrame.Histo1D("b1"); // OK, the type of "b1" is deduced at runtime
805dataFrame.Min<MyNumber_t>("myObject"); // OK, "myObject" is deduced to be of type `MyNumber_t`
806~~~
807
808Deducing types at runtime requires the just-in-time compilation of the relevant actions, which has a small runtime
809overhead, so specifying the type of the columns as template parameters to the action is good practice when performance is a goal.
810
811When deducing types at runtime, fundamental types are read as constant values, i.e. it is not possible to write to column values
812from Filters or Defines. This is typically perfectly fine and avoids certain common mistakes such as typing `x = 0` rather than `x == 0`.
813Classes and other complex types are read by non-constant references to avoid copies and to permit calls to non-const member functions.
814Note that calling non-const member functions will often not be thread-safe.
815
816### Generic actions
817`RDataFrame` strives to offer a comprehensive set of standard actions that can be performed on each event. At the same
818time, it **allows users to execute arbitrary code (i.e. a generic action) inside the event loop** through the `Foreach`
819and `ForeachSlot` actions.
820
821`Foreach(f, columnList)` takes a function `f` (lambda expression, free function, functor...) and a list of columns, and
822executes `f` on those columns for each event. The function passed must return nothing (i.e. `void`). It can be used to
823perform actions that are not already available in the interface. For example, the following snippet evaluates the root
824mean square of column "b":
825~~~{.cpp}
826// Single-thread evaluation of RMS of column "b" using Foreach
827double sumSq = 0.;
828unsigned int n = 0;
829RDataFrame d("bTree", bFilePtr);
830d.Foreach([&sumSq, &n](double b) { ++n; sumSq += b*b; }, {"b"});
831std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
832~~~
833When executing on multiple threads, users are responsible for the thread-safety of the expression passed to `Foreach`:
834each thread will execute the expression multiple times (once per entry) in an unspecified order.
835The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of `rms`; but
836this is probably too much head-scratch for such a simple operation.
837
838`ForeachSlot` can help in this situation. It is an alternative version of `Foreach` for which the function takes an
839additional parameter besides the columns it should be applied to: an `unsigned int slot` parameter, where `slot` is a
840number indicating which thread (0, 1, 2 , ..., poolSize - 1) the function is being run in. More specifically, RDataFrame
841guarantees that `ForeachSlot` will invoke the user expression with different `slot` parameters for different concurrent
842executions (there is no guarantee that a certain slot number will always correspond to a given thread id, though).
843We can take advantage of `ForeachSlot` to evaluate a thread-safe root mean square of branch "b":
844~~~{.cpp}
845// Thread-safe evaluation of RMS of branch "b" using ForeachSlot
846ROOT::EnableImplicitMT();
847const unsigned int nSlots = ROOT::GetThreadPoolSize();
848std::vector<double> sumSqs(nSlots, 0.);
849std::vector<unsigned int> ns(nSlots, 0);
850
851RDataFrame d("bTree", bFilePtr);
852d.ForeachSlot([&sumSqs, &ns](unsigned int slot, double b) { sumSqs[slot] += b*b; ns[slot] += 1; }, {"b"});
853double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
854unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
855std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
856~~~
857You see how we created one `double` variable for each thread in the pool, and later merged their results via
858`std::accumulate`.
859
860### Friend trees
861Friend trees are supported by RDataFrame.
862In order to deal with friend trees with RDataFrame, the user is required to build
863the tree and its friends and instantiate a RDataFrame with it.
864~~~{.cpp}
865TTree t([...]);
866TTree ft([...]);
867t.AddFriend(ft, "myFriend");
868
869RDataFrame d(t);
870auto f = d.Filter("myFriend.MyCol == 42");
871~~~
872
873### Reading file formats different from ROOT's
874RDataFrame can be interfaced with RDataSources. The RDataSource interface defines an API that RDataFrame can use to read arbitrary data formats.
875
876A concrete RDataSource implementation (i.e. a class that inherits from RDataSource and implements all of its pure
877methods) provides an adaptor that RDataFrame can leverage to read any kind of tabular data formats.
878RDataFrame calls into RDataSource to retrieve information about the data, retrieve (thread-local) readers or "cursors" for selected columns
879and to advance the readers to the desired data entry.
880Some predefined RDataSources are natively provided by ROOT such as the `RCsvDS` which allows to read comma separated files:
881~~~{.cpp}
882auto tdf = ROOT::RDF::MakeCsvDataFrame("MuRun2010B.csv");
883auto filteredEvents =
884 tdf.Filter("Q1 * Q2 == -1")
885 .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
886auto h = filteredEvents.Histo1D("m");
887h->Draw();
888~~~
889
890### <a name="callgraphs"></a>Call graphs (storing and reusing sets of transformations)
891**Sets of transformations can be stored as variables** and reused multiple times to create **call graphs** in which
892several paths of filtering/creation of columns are executed simultaneously; we often refer to this as "storing the
893state of the chain".
894
895This feature can be used, for example, to create a temporary column once and use it in several subsequent filters or
896actions, or to apply a strict filter to the data-set *before* executing several other transformations and actions,
897effectively reducing the amount of events processed.
898
899Let's try to make this clearer with a commented example:
900~~~{.cpp}
901// build the data-frame and specify a default column list
902RDataFrame d(treeName, filePtr, {"var1", "var2", "var3"});
903
904// apply a cut and save the state of the chain
905auto filtered = d.Filter(myBigCut);
906
907// plot branch "var1" at this point of the chain
908auto h1 = filtered.Histo1D("var1");
909
910// create a new branch "vec" with a vector extracted from a complex object (only for filtered entries)
911// and save the state of the chain
912auto newBranchFiltered = filtered.Define("vec", [](const Obj& o) { return o.getVector(); }, {"obj"});
913
914// apply a cut and fill a histogram with "vec"
915auto h2 = newBranchFiltered.Filter(cut1).Histo1D("vec");
916
917// apply a different cut and fill a new histogram
918auto h3 = newBranchFiltered.Filter(cut2).Histo1D("vec");
919
920// Inspect results
921h2->Draw(); // first access to an action result: run event-loop!
922h3->Draw("SAME"); // event loop does not need to be run again here..
923std::cout << "Entries in h1: " << h1->GetEntries() << std::endl; // ..or here
924~~~
925`RDataFrame` detects when several actions use the same filter or the same temporary column, and **only evaluates each
926filter or temporary column once per event**, regardless of how many times that result is used down the call graph.
927Objects read from each column are **built once and never copied**, for maximum efficiency.
928When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated,
929so it might be advisable to put the strictest filters first in the chain.
930
931### <a name="representgraph"></a>Printing the computation graph
932It is possible to print the computation graph from any node to obtain a dot representation either on the standard output
933or in a file.
934
935Invoking the function ROOT::RDF::SaveGraph() on any node that is not the head node, the computation graph of the branch
936the node belongs to is printed. By using the head node, the entire computation graph is printed.
937
938Following there is an example of usage:
939~~~{.cpp}
940// First, a sample computational graph is built
941ROOT::RDataFrame df("tree", "f.root");
942
943auto df2 = df.Define("x", []() { return 1; })
944 .Filter("col0 % 1 == col0")
945 .Filter([](int b1) { return b1 <2; }, {"cut1"})
946 .Define("y", []() { return 1; });
947
948auto count = df2.Count();
949
950// Prints the graph to the rd1.dot file in the current directory
951ROOT::RDF::SaveGraph(rd1, "./mydot.dot");
952// Prints the graph to standard output
953ROOT::RDF::SaveGraph(rd1);
954~~~
955
956### RDataFrame variables as function arguments and return values
957RDataFrame variables/nodes are relatively cheap to copy and it's possible to both pass them to (or move them into)
958functions and to return them from functions. However, in general each dataframe node will have a different C++ type,
959which includes all available compile-time information about what that node does. One way to cope with this complication
960is to use template functions and/or C++14 auto return types:
961~~~{.cpp}
962template <typename RDF>
963auto ApplySomeFilters(RDF df)
964{
965 return df.Filter("x > 0").Filter([](int y) { return y < 0; }, {"y"});
966}
967~~~
968
969A possibly simpler, C++11-compatible alternative is to take advantage of the fact that any dataframe node can be
970converted to the common type ROOT::RDF::RNode:
971~~~{.cpp}
972// a function that conditionally adds a Range to a RDataFrame node.
973RNode MaybeAddRange(RNode df, bool mustAddRange)
974{
975 return mustAddRange ? df.Range(1) : df;
976}
977// use as :
978ROOT::RDataFrame df(10);
979auto maybeRangedDF = MaybeAddRange(df, true);
980~~~
981
982The conversion to ROOT::RDF::RNode is cheap, but it will introduce an extra virtual call during the RDataFrame event
983loop (in most cases, the resulting performance impact should be negligible).
984
985As a final note, remember that RDataFrame actions do not return another dataframe, but a RResultPtr<T>, where T is the
986type of the result of the action.
987
988Read more on this topic [here](https://root.cern.ch/doc/master/classROOT_1_1RDF_1_1RInterface.html#a6909f04c05723de79f97a14b092318b1).
989
990<a name="reference"></a>
991*/
992// clang-format on
993
994namespace ROOT {
995
997using ColumnNamesPtr_t = std::shared_ptr<const ColumnNames_t>;
998
1000
1001////////////////////////////////////////////////////////////////////////////
1002/// \brief Build the dataframe
1003/// \param[in] treeName Name of the tree contained in the directory
1004/// \param[in] dirPtr TDirectory where the tree is stored, e.g. a TFile.
1005/// \param[in] defaultBranches Collection of default branches.
1006///
1007/// The default branches are looked at in case no branch is specified in the
1008/// booking of actions or transformations.
1009/// See RInterface for the documentation of the methods available.
1010RDataFrame::RDataFrame(std::string_view treeName, TDirectory *dirPtr, const ColumnNames_t &defaultBranches)
1011 : RInterface(std::make_shared<RDFDetail::RLoopManager>(nullptr, defaultBranches))
1012{
1013 if (!dirPtr) {
1014 auto msg = "Invalid TDirectory!";
1015 throw std::runtime_error(msg);
1016 }
1017 const std::string treeNameInt(treeName);
1018 auto tree = static_cast<TTree *>(dirPtr->Get(treeNameInt.c_str()));
1019 if (!tree) {
1020 auto msg = "Tree \"" + treeNameInt + "\" cannot be found!";
1021 throw std::runtime_error(msg);
1022 }
1023 GetProxiedPtr()->SetTree(std::shared_ptr<TTree>(tree, [](TTree *) {}));
1024}
1025
1026////////////////////////////////////////////////////////////////////////////
1027/// \brief Build the dataframe
1028/// \param[in] treeName Name of the tree contained in the directory
1029/// \param[in] filenameglob TDirectory where the tree is stored, e.g. a TFile.
1030/// \param[in] defaultBranches Collection of default branches.
1031///
1032/// The filename globbing supports the same type of expressions as TChain::Add().
1033/// The default branches are looked at in case no branch is specified in the
1034/// booking of actions or transformations.
1035/// See RInterface for the documentation of the methods available.
1036RDataFrame::RDataFrame(std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultBranches)
1037 : RInterface(std::make_shared<RDFDetail::RLoopManager>(nullptr, defaultBranches))
1038{
1039 const std::string treeNameInt(treeName);
1040 const std::string filenameglobInt(filenameglob);
1041 auto chain = std::make_shared<TChain>(treeNameInt.c_str());
1042 chain->Add(filenameglobInt.c_str());
1043 GetProxiedPtr()->SetTree(chain);
1044}
1045
1046////////////////////////////////////////////////////////////////////////////
1047/// \brief Build the dataframe
1048/// \param[in] treeName Name of the tree contained in the directory
1049/// \param[in] fileglobs Collection of file names of filename globs
1050/// \param[in] defaultBranches Collection of default branches.
1051///
1052/// The filename globbing supports the same type of expressions as TChain::Add().
1053/// The default branches are looked at in case no branch is specified in the booking of actions or transformations.
1054/// See RInterface for the documentation of the methods available.
1055RDataFrame::RDataFrame(std::string_view treeName, const std::vector<std::string> &fileglobs,
1056 const ColumnNames_t &defaultBranches)
1057 : RInterface(std::make_shared<RDFDetail::RLoopManager>(nullptr, defaultBranches))
1058{
1059 std::string treeNameInt(treeName);
1060 auto chain = std::make_shared<TChain>(treeNameInt.c_str());
1061 for (auto &f : fileglobs)
1062 chain->Add(f.c_str());
1063 GetProxiedPtr()->SetTree(chain);
1064}
1065
1066////////////////////////////////////////////////////////////////////////////
1067/// \brief Build the dataframe
1068/// \param[in] tree The tree or chain to be studied.
1069/// \param[in] defaultBranches Collection of default column names to fall back to when none is specified.
1070///
1071/// The default branches are looked at in case no branch is specified in the
1072/// booking of actions or transformations.
1073/// See RInterface for the documentation of the methods available.
1075 : RInterface(std::make_shared<RDFDetail::RLoopManager>(&tree, defaultBranches))
1076{
1077}
1078
1079//////////////////////////////////////////////////////////////////////////
1080/// \brief Build a dataframe that generates numEntries entries.
1081/// \param[in] numEntries The number of entries to generate.
1082///
1083/// An empty-source dataframe constructed with a number of entries will
1084/// generate those entries on the fly when some action is triggered,
1085/// and it will do so for all the previously-defined temporary branches.
1086/// See RInterface for the documentation of the methods available.
1088 : RInterface(std::make_shared<RDFDetail::RLoopManager>(numEntries))
1089
1090{
1091}
1092
1093//////////////////////////////////////////////////////////////////////////
1094/// \brief Build dataframe associated to datasource.
1095/// \param[in] ds The data-source object.
1096/// \param[in] defaultBranches Collection of default column names to fall back to when none is specified.
1097///
1098/// A dataframe associated to a datasource will query it to access column values.
1099/// See RInterface for the documentation of the methods available.
1100RDataFrame::RDataFrame(std::unique_ptr<ROOT::RDF::RDataSource> ds, const ColumnNames_t &defaultBranches)
1101 : RInterface(std::make_shared<RDFDetail::RLoopManager>(std::move(ds), defaultBranches))
1102{
1103}
1104
1105} // namespace ROOT
1106
1107namespace cling {
1108//////////////////////////////////////////////////////////////////////////
1109/// Print a RDataFrame at the prompt
1110std::string printValue(ROOT::RDataFrame *tdf)
1111{
1112 auto &df = *tdf->GetLoopManager();
1113 auto *tree = df.GetTree();
1114 auto defBranches = df.GetDefaultColumnNames();
1115
1116 std::ostringstream ret;
1117 if (tree) {
1118 ret << "A data frame built on top of the " << tree->GetName() << " dataset.";
1119 if (!defBranches.empty()) {
1120 if (defBranches.size() == 1)
1121 ret << "\nDefault branch: " << defBranches[0];
1122 else {
1123 ret << "\nDefault branches:\n";
1124 for (auto &&branch : defBranches) {
1125 ret << " - " << branch << "\n";
1126 }
1127 }
1128 }
1129 } else if (auto ds = tdf->fDataSource) {
1130 ret << "A data frame associated to the data source \"" << cling::printValue(ds) << "\"";
1131 } else {
1132 ret << "An empty data frame that will create " << df.GetNEmptyEntries() << " entries\n";
1133 }
1134
1135 return ret.str();
1136}
1137} // namespace cling
#define f(i)
Definition RSha256.hxx:104
unsigned long long ULong64_t
Definition RtypesCore.h:74
The head node of a RDF computation graph.
RLoopManager * GetLoopManager() const
RDataSource * fDataSource
Non-owning pointer to a data-source object. Null if no data-source. RLoopManager has ownership of the...
const std::shared_ptr< RDFDetail::RLoopManager > & GetProxiedPtr() const
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
RDataFrame(std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultBranches={})
Build the dataframe.
RDFDetail::ColumnNames_t ColumnNames_t
Describe directory structure in memory.
Definition TDirectory.h:45
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
A TTree represents a columnar dataset.
Definition TTree.h:79
std::vector< std::string > ColumnNames_t
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
std::shared_ptr< const ColumnNames_t > ColumnNamesPtr_t
Definition tree.py:1