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