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
12#include "ROOT/RDataFrame.hxx"
13#include "ROOT/RDataSource.hxx"
17#include "ROOT/RDF/Utils.hxx"
18#include <string_view>
19#include "TChain.h"
20#include "TDirectory.h"
21#include "RtypesCore.h" // for ULong64_t
22#include "TTree.h"
23
24#include <fstream> // std::ifstream
25#include <nlohmann/json.hpp> // nlohmann::json::parse
26#include <memory> // for make_shared, allocator, shared_ptr
27#include <ostream> // ostringstream
28#include <stdexcept>
29#include <string>
30#include <vector>
31
32// clang-format off
33/**
34* \class ROOT::RDataFrame
35* \ingroup dataframe
36* \brief ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree , CSV and other data formats, in C++ or Python.
37
38In addition, multi-threading and other low-level optimisations allow users to exploit all the resources available
39on their machines completely transparently.<br>
40Skip to the [class reference](#reference) or keep reading for the user guide.
41
42In a nutshell:
43~~~{.cpp}
44ROOT::EnableImplicitMT(); // Tell ROOT you want to go parallel
45ROOT::RDataFrame d("myTree", "file_*.root"); // Interface to TTree and TChain
46auto myHisto = d.Histo1D("Branch_A"); // This books the (lazy) filling of a histogram
47myHisto->Draw(); // Event loop is run here, upon first access to a result
48~~~
49
50Calculations are expressed in terms of a type-safe *functional chain of actions and transformations*, RDataFrame takes
51care of their execution. The implementation automatically puts in place several low level optimisations such as
52multi-thread parallelization and caching.
53
54\htmlonly
55<a href="https://doi.org/10.5281/zenodo.260230"><img src="https://zenodo.org/badge/DOI/10.5281/zenodo.260230.svg"
56alt="DOI"></a>
57\endhtmlonly
58
59## For the impatient user
60You can directly see RDataFrame in action in our [tutorials](https://root.cern/doc/master/group__tutorial__dataframe.html), in C++ or Python.
61
62## Table of Contents
63- [Cheat sheet](\ref cheatsheet)
64- [Introduction](\ref introduction)
65- [Crash course](\ref crash-course)
66- [Working with collections](\ref collections)
67- [Transformations: manipulating data](\ref transformations)
68- [Actions: getting results](\ref actions)
69- [Distributed execution in Python](\ref distrdf)
70- [Performance tips and parallel execution](\ref parallel-execution)
71- [More features](\ref more-features)
72 - [Systematic variations](\ref systematics)
73 - [RDataFrame objects as function arguments and return values](\ref rnode)
74 - [Storing RDataFrame objects in collections](\ref RDFCollections)
75 - [Executing callbacks every N events](\ref callbacks)
76 - [Default column lists](\ref default-branches)
77 - [Special helper columns: `rdfentry_` and `rdfslot_`](\ref helper-cols)
78 - [Just-in-time compilation: column type inference and explicit declaration of column types](\ref jitting)
79 - [User-defined custom actions](\ref generic-actions)
80 - [Dataset joins with friend trees](\ref friends)
81 - [Reading data formats other than ROOT trees](\ref other-file-formats)
82 - [Computation graphs (storing and reusing sets of transformations)](\ref callgraphs)
83 - [Visualizing the computation graph](\ref representgraph)
84 - [Activating RDataFrame execution logs](\ref rdf-logging)
85 - [Creating an RDataFrame from a dataset specification file](\ref rdf-from-spec)
86 - [Adding a progress bar](\ref progressbar)
87- [Efficient analysis in Python](\ref python)
88- <a class="el" href="classROOT_1_1RDataFrame.html#reference" onclick="javascript:toggleInherit('pub_methods_classROOT_1_1RDF_1_1RInterface')">Class reference</a>
89
90\anchor cheatsheet
91## Cheat sheet
92These are the operations which can be performed with RDataFrame.
93
94### Transformations
95Transformations are a way to manipulate the data.
96
97| **Transformation** | **Description** |
98|------------------|--------------------|
99| Alias() | Introduce an alias for a particular column name. |
100| Define() | Create 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). |
101| DefinePerSample() | Define a new column that is updated when the input sample changes, e.g. when switching tree being processed in a chain. |
102| 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`. |
103| DefineSlotEntry() | Same as DefineSlot(), but the entry number is passed in addition to the slot number. This is meant as a helper in case the expression depends on the entry number. For details about entry numbers in multi-threaded runs, see [here](\ref helper-cols). |
104| Filter() | Filter rows based on user-defined conditions. |
105| Range() | Filter rows based on entry number (single-thread only). |
106| Redefine() | Overwrite the value and/or type of an existing column. See Define() for more information. |
107| RedefineSlot() | Overwrite the value and/or type of an existing column. See DefineSlot() for more information. |
108| RedefineSlotEntry() | Overwrite the value and/or type of an existing column. See DefineSlotEntry() for more information. |
109| Vary() | Register systematic variations for an existing column. Varied results are then extracted via VariationsFor(). |
110
111
112### Actions
113Actions aggregate data into a result. Each one is described in more detail in the reference guide.
114
115In 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.
116
117Lazy actions only trigger the event loop when one of the results is accessed for the first time, making it easy to
118produce many different results in one event loop. Instant actions trigger the event loop instantly.
119
120
121| **Lazy action** | **Description** |
122|------------------|-----------------|
123| Aggregate() | Execute a user-defined accumulation operation on the processed column values. |
124| Book() | Book execution of a custom action using a user-defined helper object. |
125| Cache() | Cache column values in memory. Custom columns can be cached as well, filtered entries are not cached. Users can specify which columns to save (default is all). |
126| Count() | Return the number of events processed. Useful e.g. to get a quick count of the number of events passing a Filter. |
127| Display() | Provides a printable representation of the dataset contents. The method returns a ROOT::RDF::RDisplay() instance which can print a tabular representation of the data or return it as a string. |
128| Fill() | Fill a user-defined object with the values of the specified columns, as if by calling `Obj.Fill(col1, col2, ...)`. |
129| Graph() | Fills a TGraph with the two columns provided. If multi-threading is enabled, the order of the points may not be the one expected, it is therefore suggested to sort if before drawing. |
130| GraphAsymmErrors() | Fills a TGraphAsymmErrors. If multi-threading is enabled, the order of the points may not be the one expected, it is therefore suggested to sort if before drawing. |
131| Histo1D(), Histo2D(), Histo3D() | Fill a one-, two-, three-dimensional histogram with the processed column values. |
132| HistoND() | Fill an N-dimensional histogram with the processed column values. |
133| 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.|
134| Mean() | Return the mean of processed column values.|
135| 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.|
136| Profile1D(), Profile2D() | Fill a one- or two-dimensional profile with the column values that passed all filters. |
137| 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. |
138| Report() | Obtain 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 ROOT::RDF::RCutFlowReport instance which can be queried programmatically to get information about the effects of the individual cuts. |
139| Stats() | Return a TStatistic object filled with the input columns. |
140| StdDev() | Return the unbiased standard deviation of the processed column values. |
141| 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. |
142| Take() | Extract a column from the dataset as a collection of values, e.g. a `std::vector<float>` for a column of type `float`. |
143
144| **Instant action** | **Description** |
145|---------------------|-----------------|
146| Foreach() | Execute a user-defined function on each entry. Users are responsible for the thread-safety of this callable when executing with implicit multi-threading enabled. |
147| 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`. |
148| Snapshot() | Write the processed dataset 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 flag in the snapshot options.|
149
150
151### Queries
152
153These operations do not modify the dataframe or book computations but simply return information on the RDataFrame object.
154
155| **Operation** | **Description** |
156|---------------------|-----------------|
157| Describe() | Get useful information describing the dataframe, e.g. columns and their types. |
158| GetColumnNames() | Get the names of all the available columns of the dataset. |
159| GetColumnType() | Return the type of a given column as a string. |
160| GetColumnTypeNamesList() | Return the list of type names of columns in the dataset. |
161| GetDefinedColumnNames() | Get the names of all the defined columns. |
162| GetFilterNames() | Return the names of all filters in the computation graph. |
163| GetNRuns() | Return the number of event loops run by this RDataFrame instance so far. |
164| GetNSlots() | Return the number of processing slots that RDataFrame will use during the event loop (i.e. the concurrency level). |
165| SaveGraph() | Store the computation graph of an RDataFrame in [DOT format (graphviz)](https://en.wikipedia.org/wiki/DOT_(graph_description_language)) for easy inspection. See the [relevant section](\ref representgraph) for details. |
166
167\anchor introduction
168## Introduction
169Users define their analysis as a sequence of operations to be performed on the dataframe object; the framework
170takes care of the management of the loop over entries as well as low-level details such as I/O and parallelization.
171RDataFrame provides methods to perform most common operations required by ROOT analyses;
172at the same time, users can just as easily specify custom code that will be executed in the event loop.
173
174RDataFrame is built with a *modular* and *flexible* workflow in mind, summarised as follows:
175
1761. Construct a dataframe object by specifying a dataset. RDataFrame supports TTree as well as TChain, [CSV files](https://root.cern/doc/master/df014__CSVDataSource_8C.html), [SQLite files](https://root.cern/doc/master/df027__SQliteDependencyOverVersion_8C.html), [RNTuples](https://root.cern/doc/master/structROOT_1_1Experimental_1_1RNTuple.html), and it can be extended to custom data formats. From Python, [NumPy arrays can be imported into RDataFrame](https://root.cern/doc/master/df032__MakeNumpyDataFrame_8py.html) as well.
177
1782. Transform the dataframe by:
179
180 - [Applying filters](https://root.cern/doc/master/classROOT_1_1RDataFrame.html#transformations). This selects only specific rows of the dataset.
181
182 - [Creating custom columns](https://root.cern/doc/master/classROOT_1_1RDataFrame.html#transformations). Custom columns can, for example, contain the results of a computation that must be performed for every row of the dataset.
183
1843. [Produce results](https://root.cern/doc/master/classROOT_1_1RDataFrame.html#actions). *Actions* are used to aggregate data into results. Most actions are *lazy*, i.e. they are not executed on the spot, but registered with RDataFrame and executed only when a result is accessed for the first time.
185
186Make 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.
187
188The following table shows how analyses based on TTreeReader and TTree::Draw() translate to RDataFrame. Follow the
189[crash course](#crash-course) to discover more idiomatic and flexible ways to express analyses with RDataFrame.
190<table>
191<tr>
192 <td>
193 <b>TTreeReader</b>
194 </td>
195 <td>
196 <b>ROOT::RDataFrame</b>
197 </td>
198</tr>
199<tr>
200 <td>
201~~~{.cpp}
202TTreeReader reader("myTree", file);
203TTreeReaderValue<A_t> a(reader, "A");
204TTreeReaderValue<B_t> b(reader, "B");
205TTreeReaderValue<C_t> c(reader, "C");
206while(reader.Next()) {
207 if(IsGoodEvent(*a, *b, *c))
208 DoStuff(*a, *b, *c);
209}
210~~~
211 </td>
212 <td>
213~~~{.cpp}
214ROOT::RDataFrame d("myTree", file, {"A", "B", "C"});
215d.Filter(IsGoodEvent).Foreach(DoStuff);
216~~~
217 </td>
218</tr>
219<tr>
220 <td>
221 <b>TTree::Draw</b>
222 </td>
223 <td>
224 <b>ROOT::RDataFrame</b>
225 </td>
226</tr>
227<tr>
228 <td>
229~~~{.cpp}
230auto *tree = file->Get<TTree>("myTree");
231tree->Draw("x", "y > 2");
232~~~
233 </td>
234 <td>
235~~~{.cpp}
236ROOT::RDataFrame df("myTree", file);
237auto h = df.Filter("y > 2").Histo1D("x");
238h->Draw()
239~~~
240 </td>
241</tr>
242<tr>
243 <td>
244~~~{.cpp}
245tree->Draw("jet_eta", "weight*(event == 1)");
246~~~
247 </td>
248 <td>
249~~~{.cpp}
250df.Filter("event == 1").Histo1D("jet_eta", "weight");
251// or the fully compiled version:
252df.Filter([] (ULong64_t e) { return e == 1; }, {"event"}).Histo1D<RVec<float>>("jet_eta", "weight");
253~~~
254 </td>
255</tr>
256<tr>
257 <td>
258~~~{cpp}
259// object selection: for each event, fill histogram with array of selected pts
260tree->Draw('Muon_pt', 'Muon_pt > 100')
261~~~
262 </td>
263 <td>
264~~~{cpp}
265// with RDF, arrays are read as ROOT::VecOps::RVec objects
266df.Define("good_pt", "Muon_pt[Muon_pt > 100]").Histo1D("good_pt")
267~~~
268 </td>
269</tr>
270</table>
271
272\anchor crash-course
273## Crash course
274All snippets of code presented in the crash course can be executed in the ROOT interpreter. Simply precede them with
275~~~{.cpp}
276using namespace ROOT; // RDataFrame's namespace
277~~~
278which is omitted for brevity. The terms "column" and "branch" are used interchangeably.
279
280### Creating an RDataFrame
281RDataFrame's constructor is where the user specifies the dataset and, optionally, a default set of columns that
282operations should work with. Here are the most common methods to construct an RDataFrame object:
283~~~{.cpp}
284// single file -- all constructors are equivalent
285TFile *f = TFile::Open("file.root");
286TTree *t = f.Get<TTree>("treeName");
287
288RDataFrame d1("treeName", "file.root");
289RDataFrame d2("treeName", f); // same as TTreeReader
290RDataFrame d3(*t);
291
292// multiple files -- all constructors are equivalent
293TChain chain("myTree");
294chain.Add("file1.root");
295chain.Add("file2.root");
296
297RDataFrame d4("myTree", {"file1.root", "file2.root"});
298std::vector<std::string> files = {"file1.root", "file2.root"};
299RDataFrame d5("myTree", files);
300RDataFrame d6("myTree", "file*.root"); // the glob is passed as-is to TChain's constructor
301RDataFrame d7(chain);
302~~~
303Additionally, users can construct an RDataFrame with no data source by passing an integer number. This is the number of rows that
304will be generated by this RDataFrame.
305~~~{.cpp}
306RDataFrame d(10); // a RDF with 10 entries (and no columns/branches, for now)
307d.Foreach([] { static int i = 0; std::cout << i++ << std::endl; }); // silly example usage: count to ten
308~~~
309This is useful to generate simple datasets 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](https://pythia.org/) events and write them to disk in parallel (with the Snapshot action).
310
311For data sources other than TTrees and TChains, RDataFrame objects are constructed using ad-hoc factory functions (see e.g. FromCSV(), FromSqlite(), FromArrow()):
312
313~~~{.cpp}
314auto df = ROOT::RDF::FromCSV("input.csv");
315// use df as usual
316~~~
317
318### Filling a histogram
319Let's now tackle a very common task, filling a histogram:
320~~~{.cpp}
321// Fill a TH1D with the "MET" branch
322RDataFrame d("myTree", "file.root");
323auto h = d.Histo1D("MET");
324h->Draw();
325~~~
326The first line creates an RDataFrame associated to the TTree "myTree". This tree has a branch named "MET".
327
328Histo1D() is an *action*; it returns a smart pointer (a ROOT::RDF::RResultPtr, to be precise) to a TH1D histogram filled
329with the `MET` of all events. If the quantity stored in the column is a collection (e.g. a vector or an array), the
330histogram is filled with all vector elements for each event.
331
332You can use the objects returned by actions as if they were pointers to the desired results. There are many other
333possible [actions](\ref cheatsheet), and all their results are wrapped in smart pointers; we'll see why in a minute.
334
335### Applying a filter
336Let'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:
337~~~{.cpp}
338RDataFrame d("myTree", "file.root");
339auto c = d.Filter("MET > 4.").Count(); // computations booked, not run
340std::cout << *c << std::endl; // computations run here, upon first access to the result
341~~~
342The filter string (which must contain a valid C++ expression) is applied to the specified columns for each event;
343the name and types of the columns are inferred automatically. The string expression is required to return a `bool`
344which signals whether the event passes the filter (`true`) or not (`false`).
345
346You can think of your data as "flowing" through the chain of calls, being transformed, filtered and finally used to
347perform actions. Multiple Filter() calls can be chained one after another.
348
349Using string filters is nice for simple things, but they are limited to specifying the equivalent of a single return
350statement or the body of a lambda, so it's cumbersome to use strings with more complex filters. They also add a small
351runtime overhead, as ROOT needs to just-in-time compile the string into C++ code. When more freedom is required or
352runtime performance is very important, a C++ callable can be specified instead (a lambda in the following snippet,
353but it can be any kind of function or even a functor class), together with a list of column names.
354This snippet is analogous to the one above:
355~~~{.cpp}
356RDataFrame d("myTree", "file.root");
357auto metCut = [](double x) { return x > 4.; }; // a C++11 lambda function checking "x > 4"
358auto c = d.Filter(metCut, {"MET"}).Count();
359std::cout << *c << std::endl;
360~~~
361
362An example of a more complex filter expressed as a string containing C++ code is shown below
363
364~~~{.cpp}
365RDataFrame d("myTree", "file.root");
366auto df = d.Define("p", "std::array<double, 4> p{px, py, pz}; return p;")
367 .Filter("double p2 = 0.0; for (auto&& x : p) p2 += x*x; return sqrt(p2) < 10.0;");
368~~~
369
370The code snippet above defines a column `p` that is a fixed-size array using the component column names and then
371filters on its magnitude by looping over its elements. It must be noted that the usage of strings to define columns
372like the one above is currently the only possibility when using PyROOT. When writing expressions as such, only constants
373and data coming from other columns in the dataset can be involved in the code passed as a string. Local variables and
374functions cannot be used, since the interpreter will not know how to find them. When capturing local state is necessary,
375it must first be declared to the ROOT C++ interpreter.
376
377More information on filters and how to use them to automatically generate cutflow reports can be found [below](#Filters).
378
379### Defining custom columns
380Let's now consider the case in which "myTree" contains two quantities "x" and "y", but our analysis relies on a derived
381quantity `z = sqrt(x*x + y*y)`. Using the Define() transformation, we can create a new column in the dataset containing
382the variable "z":
383~~~{.cpp}
384RDataFrame d("myTree", "file.root");
385auto sqrtSum = [](double x, double y) { return sqrt(x*x + y*y); };
386auto zMean = d.Define("z", sqrtSum, {"x","y"}).Mean("z");
387std::cout << *zMean << std::endl;
388~~~
389Define() creates the variable "z" by applying `sqrtSum` to "x" and "y". Later in the chain of calls we refer to
390variables created with Define() as if they were actual tree branches/columns, but they are evaluated on demand, at most
391once per event. As with filters, Define() calls can be chained with other transformations to create multiple custom
392columns. Define() and Filter() transformations can be concatenated and intermixed at will.
393
394As with filters, it is possible to specify new columns as string expressions. This snippet is analogous to the one above:
395~~~{.cpp}
396RDataFrame d("myTree", "file.root");
397auto zMean = d.Define("z", "sqrt(x*x + y*y)").Mean("z");
398std::cout << *zMean << std::endl;
399~~~
400
401Again the names of the columns used in the expression and their types are inferred automatically. The string must be
402valid C++ and it is just-in-time compiled. The process has a small runtime overhead and like with filters it is currently the only possible approach when using PyROOT.
403
404Previously, when showing the different ways an RDataFrame can be created, we showed a constructor that takes a
405number of entries as a parameter. In the following example we show how to combine such an "empty" RDataFrame with Define()
406transformations to create a dataset on the fly. We then save the generated data on disk using the Snapshot() action.
407~~~{.cpp}
408RDataFrame d(100); // an RDF that will generate 100 entries (currently empty)
409int x = -1;
410auto d_with_columns = d.Define("x", [&x] { return ++x; })
411 .Define("xx", [&x] { return x*x; });
412d_with_columns.Snapshot("myNewTree", "newfile.root");
413~~~
414This example is slightly more advanced than what we have seen so far. First, it makes use of lambda captures (a
415simple way to make external variables available inside the body of C++ lambdas) to act on the same variable `x` from
416both Define() transformations. Second, we have *stored* the transformed dataframe in a variable. This is always
417possible, since at each point of the transformation chain users can store the status of the dataframe for further use (more
418on this [below](#callgraphs)).
419
420You can read more about defining new columns [here](#custom-columns).
421
422\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."
423
424
425### Running on a range of entries
426It is sometimes necessary to limit the processing of the dataset to a range of entries. For this reason, the RDataFrame
427offers the concept of ranges as a node of the RDataFrame chain of transformations; this means that filters, columns and
428actions can be concatenated to and intermixed with Range()s. If a range is specified after a filter, the range will act
429exclusively on the entries passing the filter -- it will not even count the other entries! The same goes for a Range()
430hanging from another Range(). Here are some commented examples:
431~~~{.cpp}
432RDataFrame d("myTree", "file.root");
433// Here we store a dataframe that loops over only the first 30 entries in a variable
434auto d30 = d.Range(30);
435// This is how you pick all entries from 15 onwards
436auto d15on = d.Range(15, 0);
437// We can specify a stride too, in this case we pick an event every 3
438auto d15each3 = d.Range(0, 15, 3);
439~~~
440Note that ranges are not available when multi-threading is enabled. More information on ranges is available
441[here](#ranges).
442
443### Executing multiple actions in the same event loop
444As a final example let us apply two different cuts on branch "MET" and fill two different histograms with the "pt_v" of
445the filtered events.
446By now, you should be able to easily understand what is happening:
447~~~{.cpp}
448RDataFrame d("treeName", "file.root");
449auto h1 = d.Filter("MET > 10").Histo1D("pt_v");
450auto h2 = d.Histo1D("pt_v");
451h1->Draw(); // event loop is run once here
452h2->Draw("SAME"); // no need to run the event loop again
453~~~
454RDataFrame executes all above actions by **running the event-loop only once**. The trick is that actions are not
455executed at the moment they are called, but they are **lazy**, i.e. delayed until the moment one of their results is
456accessed through the smart pointer. At that time, the event loop is triggered and *all* results are produced
457simultaneously.
458
459It is therefore good practice to declare all your transformations and actions *before* accessing their results, allowing
460RDataFrame to run the loop once and produce all results in one go.
461
462### Going parallel
463Let's say we would like to run the previous examples in parallel on several cores, dividing events fairly between cores.
464The only modification required to the snippets would be the addition of this line *before* constructing the main
465dataframe object:
466~~~{.cpp}
467ROOT::EnableImplicitMT();
468~~~
469Simple as that. More details are given [below](#parallel-execution).
470
471\anchor collections
472## Working with collections and object selections
473
474RDataFrame reads collections as the special type [ROOT::RVec](https://root.cern/doc/master/classROOT_1_1VecOps_1_1RVec.html): for example, a column containing an array of floating point numbers can be read as a ROOT::RVecF. C-style arrays (with variable or static size), STL vectors and most other collection types can be read this way.
475
476RVec 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 vectorised fashion, similarly to Python's NumPy arrays.
477
478For example, to fill a histogram with the "pt" of selected particles for each event, Define() can be used to create a column that contains the desired array elements as follows:
479
480~~~{.cpp}
481// h is filled with all the elements of `good_pts`, for each event
482auto h = df.Define("good_pts", [](const ROOT::RVecF &pt) { return pt[pt > 0]; })
483 .Histo1D("good_pts");
484~~~
485
486And in Python:
487
488~~~{.py}
489h = df.Define("good_pts", "pt[pt > 0]").Histo1D("good_pts")
490~~~
491
492Learn more at ROOT::VecOps::RVec.
493
494\anchor transformations
495## Transformations: manipulating data
496\anchor Filters
497### Filters
498A filter is created through a call to `Filter(f, columnList)` or `Filter(filterString)`. In the first overload, `f` can
499be a function, a lambda expression, a functor class, or any other callable object. It must return a `bool` signalling
500whether the event has passed the selection (`true`) or not (`false`). It should perform "read-only" operations on the
501columns, and should not have side-effects (e.g. modification of an external or static variable) to ensure correctness
502when implicit multi-threading is active. The second overload takes a string with a valid C++ expression in which column
503names are used as variable names (e.g. `Filter("x[0] + x[1] > 0")`). This is a convenience feature that comes with a
504certain runtime overhead: C++ code has to be generated on the fly from this expression before using it in the event
505loop. See the paragraph about "Just-in-time compilation" below for more information.
506
507RDataFrame only evaluates filters when necessary: if multiple filters are chained one after another, they are executed
508in order and the first one returning `false` causes the event to be discarded and triggers the processing of the next
509entry. If multiple actions or transformations depend on the same filter, that filter is not executed multiple times for
510each entry: after the first access it simply serves a cached result.
511
512\anchor named-filters-and-cutflow-reports
513#### Named filters and cutflow reports
514An optional string parameter `name` can be passed to the Filter() method to create a **named filter**. Named filters
515work as usual, but also keep track of how many entries they accept and reject.
516
517Statistics are retrieved through a call to the Report() method:
518
519- when Report() is called on the main RDataFrame object, it returns a ROOT::RDF::RResultPtr<RCutFlowReport> relative to all
520named filters declared up to that point
521- when called on a specific node (e.g. the result of a Define() or Filter()), it returns a ROOT::RDF::RResultPtr<RCutFlowReport>
522relative all named filters in the section of the chain between the main RDataFrame and that node (included).
523
524Stats are stored in the same order as named filters have been added to the graph, and *refer to the latest event-loop*
525that has been run using the relevant RDataFrame.
526
527\anchor ranges
528### Ranges
529When RDataFrame is not being used in a multi-thread environment (i.e. no call to EnableImplicitMT() was made),
530Range() transformations are available. These act very much like filters but instead of basing their decision on
531a filter expression, they rely on `begin`,`end` and `stride` parameters.
532
533- `begin`: initial entry number considered for this range.
534- `end`: final entry number (excluded) considered for this range. 0 means that the range goes until the end of the dataset.
535- `stride`: process one entry of the [begin, end) range every `stride` entries. Must be strictly greater than 0.
536
537The actual number of entries processed downstream of a Range() node will be `(end - begin)/stride` (or less if less
538entries than that are available).
539
540Note that ranges act "locally", not based on the global entry count: `Range(10,50)` means "skip the first 10 entries
541*that reach this node*, let the next 40 entries pass, then stop processing". If a range node hangs from a filter node,
542and the range has a `begin` parameter of 10, that means the range will skip the first 10 entries *that pass the
543preceding filter*.
544
545Ranges allow "early quitting": if all branches of execution of a functional graph reached their `end` value of
546processed entries, the event-loop is immediately interrupted. This is useful for debugging and quick data explorations.
547
548\anchor custom-columns
549### Custom columns
550Custom columns are created by invoking `Define(name, f, columnList)`. As usual, `f` can be any callable object
551(function, lambda expression, functor class...); it takes the values of the columns listed in `columnList` (a list of
552strings) as parameters, in the same order as they are listed in `columnList`. `f` must return the value that will be
553assigned to the temporary column.
554
555A new variable is created called `name`, accessible as if it was contained in the dataset from subsequent
556transformations/actions.
557
558Use cases include:
559- caching the results of complex calculations for easy and efficient multiple access
560- extraction of quantities of interest from complex objects
561- branch aliasing, i.e. changing the name of a branch
562
563An exception is thrown if the `name` of the new column/branch is already in use for another branch in the TTree.
564
565It is also possible to specify the quantity to be stored in the new temporary column as a C++ expression with the method
566`Define(name, expression)`. For example this invocation
567
568~~~{.cpp}
569df.Define("pt", "sqrt(px*px + py*py)");
570~~~
571
572will create a new column called "pt" the value of which is calculated starting from the columns px and py. The system
573builds a just-in-time compiled function starting from the expression after having deduced the list of necessary branches
574from the names of the variables specified by the user.
575
576#### Custom columns as function of slot and entry number
577
578It is possible to create custom columns also as a function of the processing slot and entry numbers. The methods that can
579be invoked are:
580- `DefineSlot(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, T1, T2, ...)`: the
581first parameter is the slot number which ranges from 0 to ROOT::GetThreadPoolSize() - 1.
582- `DefineSlotEntry(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, ULong64_t,
583T1, T2, ...)`: the first parameter is the slot number while the second one the number of the entry being processed.
584
585\anchor actions
586## Actions: getting results
587### Instant and lazy actions
588Actions can be **instant** or **lazy**. Instant actions are executed as soon as they are called, while lazy actions are
589executed whenever the object they return is accessed for the first time. As a rule of thumb, actions with a return value
590are lazy, the others are instant.
591
592### Return type of a lazy action
593
594When a lazy action is called, it returns a \link ROOT::RDF::RResultPtr ROOT::RDF::RResultPtr<T>\endlink, where T is the
595type of the result of the action. The final result will be stored in the `RResultPtr` and can be retrieved by
596dereferencing it or via its `GetValue` method.
597
598### Actions that return collections
599
600If the type of the return value of an action is a collection, e.g. `std::vector<int>`, you can iterate its elements
601directly through the wrapping `RResultPtr`:
602
603~~~{.cpp}
604ROOT::RDataFrame df{5};
605auto df1 = df.Define("x", []{ return 42; });
606for (const auto &el: df1.Take<int>("x")){
607 std::cout << "Element: " << el << "\n";
608}
609~~~
610
611~~~{.py}
612df = ROOT.RDataFrame(5).Define("x", "42")
613for el in df.Take[int]("x"):
614 print(f"Element: {el}")
615~~~
616
617### Actions and readers
618
619An action that needs values for its computations will request it from a reader, e.g. a column created via `Define` or
620available from the input dataset. The action will request values from each column of the list of input columns (either
621inferred or specified by the user), in order. For example:
622
623~~~{.cpp}
624ROOT::RDataFrame df{1};
625auto df1 = df.Define("x", []{ return 11; });
626auto df2 = df1.Define("y", []{ return 22; });
627auto graph = df2.Graph<int, int>("x","y");
628~~~
629
630The `Graph` action is going to request first the value from column "x", then that of column "y". Specifically, the order
631of execution of the operations of nodes in this branch of the computation graph is guaranteed to be top to bottom.
632
633\anchor distrdf
634## Distributed execution
635
636RDataFrame applications can be executed in parallel through distributed computing frameworks on a set of remote machines
637thanks to the Python package `ROOT.RDF.Experimental.Distributed`. This experimental, **Python-only** package allows to scale the
638optimized performance RDataFrame can achieve on a single machine to multiple nodes at the same time. It is designed so
639that different backends can be easily plugged in, currently supporting [Apache Spark](http://spark.apache.org/) and
640[Dask](https://dask.org/). To make use of distributed RDataFrame, you only need to switch `ROOT.RDataFrame` with
641the backend-specific `RDataFrame` of your choice, for example:
642
643~~~{.py}
644import ROOT
645
646# Point RDataFrame calls to the Spark specific RDataFrame
647RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
648
649# It still accepts the same constructor arguments as traditional RDataFrame
650df = RDataFrame("mytree", "myfile.root")
651
652# Continue the application with the traditional RDataFrame API
653sum = df.Filter("x > 10").Sum("y")
654h = df.Histo1D(("name", "title", 10, 0, 10), "x")
655
656print(sum.GetValue())
657h.Draw()
658~~~
659
660The main goal of this package is to support running any RDataFrame application distributedly. Nonetheless, not all
661parts of the RDataFrame API currently work with this package. The subset that is currently available is:
662- AsNumpy
663- Count
664- Define
665- DefinePerSample
666- Filter
667- Graph
668- Histo[1,2,3]D
669- HistoND
670- Max
671- Mean
672- Min
673- Profile[1,2,3]D
674- Redefine
675- Snapshot
676- Stats
677- StdDev
678- Sum
679- Systematic variations: Vary and [VariationsFor](\ref ROOT::RDF::Experimental::VariationsFor).
680- Parallel submission of distributed graphs: [RunGraphs](\ref ROOT::RDF::RunGraphs).
681- Information about the dataframe: GetColumnNames.
682
683with support for more operations coming in the future. Data sources other than TTree and TChain (e.g. CSV, RNTuple) are
684currently not supported.
685
686\note The distributed RDataFrame module requires at least Python version 3.8.
687
688### Connecting to a Spark cluster
689
690In order to distribute the RDataFrame workload, you can connect to a Spark cluster you have access to through the
691official [Spark API](https://spark.apache.org/docs/latest/rdd-programming-guide.html#initializing-spark), then hook the
692connection instance to the distributed `RDataFrame` object like so:
693
694~~~{.py}
695import pyspark
696import ROOT
697
698# Create a SparkContext object with the right configuration for your Spark cluster
699conf = SparkConf().setAppName(appName).setMaster(master)
700sc = SparkContext(conf=conf)
701
702# Point RDataFrame calls to the Spark specific RDataFrame
703RDataFrame = ROOT.RDF.Experimental.Distributed.Spark.RDataFrame
704
705# The Spark RDataFrame constructor accepts an optional "sparkcontext" parameter
706# and it will distribute the application to the connected cluster
707df = RDataFrame("mytree", "myfile.root", sparkcontext = sc)
708~~~
709
710If an instance of [SparkContext](https://spark.apache.org/docs/latest/api/python/reference/api/pyspark.SparkContext.html)
711is not provided, the default behaviour is to create one in the background for you.
712
713### Connecting to a Dask cluster
714
715Similarly, you can connect to a Dask cluster by creating your own connection object which internally operates with one
716of the cluster schedulers supported by Dask (more information in the
717[Dask distributed docs](http://distributed.dask.org/en/stable/)):
718
719~~~{.py}
720import ROOT
721from dask.distributed import Client
722
723# Point RDataFrame calls to the Dask specific RDataFrame
724RDataFrame = ROOT.RDF.Experimental.Distributed.Dask.RDataFrame
725
726# In a Python script the Dask client needs to be initalized in a context
727# Jupyter notebooks / Python session don't need this
728if __name__ == "__main__":
729 # With an already setup cluster that exposes a Dask scheduler endpoint
730 client = Client("dask_scheduler.domain.com:8786")
731
732 # The Dask RDataFrame constructor accepts the Dask Client object as an optional argument
733 df = RDataFrame("mytree","myfile.root", daskclient=client)
734 # Proceed as usual
735 df.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
736~~~
737
738If an instance of [distributed.Client](http://distributed.dask.org/en/stable/api.html#distributed.Client) is not
739provided to the RDataFrame object, it will be created for you and it will run the computations in the local machine
740using all cores available.
741
742### Choosing the number of distributed tasks
743
744A distributed RDataFrame has internal logic to define in how many chunks the input dataset will be split before sending
745tasks to the distributed backend. Each task reads and processes one of said chunks. The logic is backend-dependent, but
746generically tries to infer how many cores are available in the cluster through the connection object. The number of
747tasks will be equal to the inferred number of cores. There are cases where the connection object of the chosen backend
748doesn't have information about the actual resources of the cluster. An example of this is when using Dask to connect to
749a batch system. The client object created at the beginning of the application does not automatically know how many cores
750will be available during distributed execution, since the jobs are submitted to the batch system after the creation of
751the connection. In such cases, the logic is to default to process the whole dataset in 2 tasks.
752
753The number of tasks submitted for distributed execution can be also set programmatically, by providing the optional
754keyword argument `npartitions` when creating the RDataFrame object. This parameter is accepted irrespectively of the
755backend used:
756
757~~~{.py}
758import ROOT
759
760# Define correct imports and access the distributed RDataFrame appropriate for the
761# backend used in the analysis
762RDataFrame = ROOT.RDF.Experimental.Distributed.[BACKEND].RDataFrame
763
764if __name__ == "__main__":
765 # The `npartitions` optional argument tells the RDataFrame how many tasks are desired
766 df = RDataFrame("mytree","myfile.root", npartitions=NPARTITIONS)
767 # Proceed as usual
768 df.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
769~~~
770
771Note that when processing a TTree or TChain dataset, the `npartitions` value should not exceed the number of clusters in
772the dataset. The number of clusters in a TTree can be retrieved by typing `rootls -lt myfile.root` at a command line.
773
774### Distributed Snapshot
775
776The Snapshot operation behaves slightly differently when executed distributedly. First off, it requires the path
777supplied to the Snapshot call to be accessible from any worker of the cluster and from the client machine (in general
778it should be provided as an absolute path). Another important difference is that `n` separate files will be produced,
779where `n` is the number of dataset partitions. As with local RDataFrame, the result of a Snapshot on a distributed
780RDataFrame is another distributed RDataFrame on which we can define a new computation graph and run more distributed
781computations.
782
783### Distributed RunGraphs
784
785Submitting multiple distributed RDataFrame executions is supported through the RunGraphs function. Similarly to its
786local counterpart, the function expects an iterable of objects representing an RDataFrame action. Each action will be
787triggered concurrently to send multiple computation graphs to a distributed cluster at the same time:
788
789~~~{.py}
790import ROOT
791RDataFrame = ROOT.RDF.Experimental.Distributed.Dask.RDataFrame
792RunGraphs = ROOT.RDF.Experimental.Distributed.RunGraphs
793
794# Create 3 different dataframes and book an histogram on each one
795histoproxies = [
796 RDataFrame(100)
797 .Define("x", "rdfentry_")
798 .Histo1D(("name", "title", 10, 0, 100), "x")
799 for _ in range(4)
800]
801
802# Execute the 3 computation graphs
803RunGraphs(histoproxies)
804# Retrieve all the histograms in one go
805histos = [histoproxy.GetValue() for histoproxy in histoproxies]
806~~~
807
808Every distributed backend supports this feature and graphs belonging to different backends can be still triggered with
809a single call to RunGraphs (e.g. it is possible to send a Spark job and a Dask job at the same time).
810
811### Histogram models in distributed mode
812
813When calling a Histo*D operation in distributed mode, remember to pass to the function the model of the histogram to be
814computed, e.g. the axis range and the number of bins:
815
816~~~{.py}
817import ROOT
818RDataFrame = ROOT.RDF.Experimental.Distributed.[BACKEND].RDataFrame
819
820if __name__ == "__main__":
821 df = RDataFrame("mytree","myfile.root").Define("x","someoperation")
822 # The model can be passed either as a tuple with the arguments in the correct order
823 df.Histo1D(("name", "title", 10, 0, 10), "x")
824 # Or by creating the specific struct
825 model = ROOT.RDF.TH1DModel("name", "title", 10, 0, 10)
826 df.Histo1D(model, "x")
827~~~
828
829Without this, two partial histograms resulting from two distributed tasks would have incompatible binning, thus leading
830to errors when merging them. Failing to pass a histogram model will raise an error on the client side, before starting
831the distributed execution.
832
833### Live visualization in distributed mode with dask
834
835The live visualization feature allows real-time data representation of plots generated during the execution
836of a distributed RDataFrame application.
837It enables visualizing intermediate results as they are computed across multiple nodes of a Dask cluster
838by creating a canvas and continuously updating it as partial results become available.
839
840The LiveVisualize() function can be imported from the Python package **ROOT.RDF.Experimental.Distributed**:
841
842~~~{.py}
843import ROOT
844
845LiveVisualize = ROOT.RDF.Experimental.Distributed.LiveVisualize
846~~~
847
848The function takes drawable objects (e.g. histograms) and optional callback functions as argument, it accepts 4 different input formats:
849
850- Passing a list or tuple of drawables:
851You can pass a list or tuple containing the plots you want to visualize. For example:
852
853~~~{.py}
854LiveVisualize([h_gaus, h_exp, h_random])
855~~~
856
857- Passing a list or tuple of drawables with a global callback function:
858You can also include a global callback function that will be applied to all plots. For example:
859
860~~~{.py}
861def set_fill_color(hist):
862 hist.SetFillColor(ROOT.kBlue)
863
864LiveVisualize([h_gaus, h_exp, h_random], set_fill_color)
865~~~
866
867- Passing a Dictionary of drawables and callback functions:
868For more control, you can create a dictionary where keys are plots and values are corresponding (optional) callback functions. For example:
869
870~~~{.py}
871plot_callback_dict = {
872 graph: set_marker,
873 h_exp: fit_exp,
874 tprofile_2d: None
875}
876
877LiveVisualize(plot_callback_dict)
878~~~
879
880- Passing a Dictionary of drawables and callback functions with a global callback function:
881You can also combine a dictionary of plots and callbacks with a global callback function:
882
883~~~{.py}
884LiveVisualize(plot_callback_dict, write_to_tfile)
885~~~
886
887\note The allowed operations to pass to LiveVisualize are:
888 - Histo1D(), Histo2D(), Histo3D()
889 - Graph()
890 - Profile1D(), Profile2D()
891
892\warning The Live Visualization feature is only supported for the Dask backend.
893
894\anchor parallel-execution
895## Performance tips and parallel execution
896As pointed out before in this document, RDataFrame can transparently perform multi-threaded event loops to speed up
897the execution of its actions. Users have to call ROOT::EnableImplicitMT() *before* constructing the RDataFrame
898object to indicate that it should take advantage of a pool of worker threads. **Each worker thread processes a distinct
899subset of entries**, and their partial results are merged before returning the final values to the user.
900There are no guarantees on the order in which threads will process the batches of entries.
901In particular, note that this means that, for multi-thread event loops, there is no
902guarantee on the order in which Snapshot() will _write_ entries: they could be scrambled with respect to the input dataset. The values of the special `rdfentry_` column will also not correspond to the entry numbers in the input dataset (e.g. TChain) in multi-thread runs.
903
904\warning By default, RDataFrame will use as many threads as the hardware supports, using up **all** the resources on
905a machine. This might be undesirable on shared computing resources such as a batch cluster. Therefore, when running on shared computing resources, use
906~~~{.cpp}
907ROOT::EnableImplicitMT(i)
908~~~
909replacing `i` with the number of CPUs/slots that were allocated for this job.
910
911### Thread-safety of user-defined expressions
912RDataFrame operations such as Histo1D() or Snapshot() are guaranteed to work correctly in multi-thread event loops.
913User-defined expressions, such as strings or lambdas passed to Filter(), Define(), Foreach(), Reduce() or Aggregate()
914will have to be thread-safe, i.e. it should be possible to call them concurrently from different threads.
915
916Note that simple Filter() and Define() transformations will inherently satisfy this requirement: Filter() / Define()
917expressions will often be *pure* in the functional programming sense (no side-effects, no dependency on external state),
918which eliminates all risks of race conditions.
919
920In order to facilitate writing of thread-safe operations, some RDataFrame features such as Foreach(), Define() or \link ROOT::RDF::RResultPtr::OnPartialResult OnPartialResult()\endlink
921offer thread-aware counterparts (ForeachSlot(), DefineSlot(), \link ROOT::RDF::RResultPtr::OnPartialResultSlot OnPartialResultSlot()\endlink): their only difference is that they
922will pass an extra `slot` argument (an unsigned integer) to the user-defined expression. When calling user-defined code
923concurrently, RDataFrame guarantees that different threads will employ different values of the `slot` parameter,
924where `slot` will be a number between 0 and `GetNSlots() - 1`.
925In other words, within a slot, computation runs sequentially and events are processed sequentially.
926Note that the same slot might be associated to different threads over the course of a single event loop, but two threads
927will never receive the same slot at the same time.
928This extra parameter might facilitate writing safe parallel code by having each thread write/modify a different
929processing slot, e.g. a different element of a list. See [here](#generic-actions) for an example usage of ForeachSlot().
930
931### Parallel execution of multiple RDataFrame event loops
932A complex analysis may require multiple separate RDataFrame computation graphs to produce all desired results. This poses the challenge that the
933event loops of each computation graph can be parallelized, but the different loops run sequentially, one after the other.
934On many-core architectures it might be desirable to run different event loops concurrently to improve resource usage.
935ROOT::RDF::RunGraphs() allows running multiple RDataFrame event loops concurrently:
936~~~{.cpp}
937ROOT::EnableImplicitMT();
938ROOT::RDataFrame df1("tree1", "f1.root");
939ROOT::RDataFrame df2("tree2", "f2.root");
940auto histo1 = df1.Histo1D("x");
941auto histo2 = df2.Histo1D("y");
942
943// just accessing result pointers, the event loops of separate RDataFrames run one after the other
944histo1->Draw(); // runs first multi-thread event loop
945histo2->Draw(); // runs second multi-thread event loop
946
947// alternatively, with ROOT::RDF::RunGraphs, event loops for separate computation graphs can run concurrently
948ROOT::RDF::RunGraphs({histo1, histo2});
949histo1->Draw(); // results can then be used as usual
950~~~
951
952### Performance considerations
953
954To obtain the maximum performance out of RDataFrame, make sure to avoid just-in-time compiled versions of transformations and actions if at all possible.
955For instance, `Filter("x > 0")` requires just-in-time compilation of the corresponding C++ logic, while the equivalent `Filter([](float x) { return x > 0.; }, {"x"})` does not.
956Similarly, `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
957should be preferred for performance-critical applications.
958
959Python applications cannot easily specify template parameters or pass C++ callables to RDataFrame.
960See [Efficient analysis in Python](#python) for possible ways to speed up hot paths in this case.
961
962Just-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*
963before 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.
964
965Also 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). See [Activating RDataFrame execution logs](\ref rdf-logging).
966
967### Memory usage
968
969There 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. Reducing the number of threads or using coarser binning will reduce the memory usage.
970
971Secondly, 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:
972
973~~~{.cpp}
974// assuming df1 and df2 are separate computation graphs, do:
975auto h1 = df1.Histo1D("x");
976auto h2 = df2.Histo1D("y");
977h1->Draw(); // we just-in-time compile everything needed by df1 and df2 here
978h2->Draw("SAME");
979
980// do not:
981auto h1 = df1.Histo1D("x");
982h1->Draw(); // we just-in-time compile here
983auto h2 = df2.Histo1D("y");
984h2->Draw("SAME"); // we just-in-time compile again here, as the second Histo1D call is new
985~~~
986
987\anchor more-features
988## More features
989Here is a list of the most important features that have been omitted in the "Crash course" for brevity.
990You don't need to read all these to start using RDataFrame, but they are useful to save typing time and runtime.
991
992\anchor systematics
993### Systematic variations
994
995Starting from ROOT v6.26, RDataFrame provides a flexible syntax to define systematic variations.
996This is done in two steps: a) register variations for one or more existing columns using Vary() and b) extract variations
997of normal RDataFrame results using \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()". In between these steps, no other change
998to the analysis code is required: the presence of systematic variations for certain columns is automatically propagated
999through filters, defines and actions, and RDataFrame will take these dependencies into account when producing varied
1000results. \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" is included in header `ROOT/RDFHelpers.hxx`. The compiled C++ programs must include this header
1001explicitly, this is not required for ROOT macros.
1002
1003An example usage of Vary() and \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" in C++:
1004
1005~~~{.cpp}
1006auto nominal_hx =
1007 df.Vary("pt", "ROOT::RVecD{pt*0.9f, pt*1.1f}", {"down", "up"})
1008 .Filter("pt > pt_cut")
1009 .Define("x", someFunc, {"pt"})
1010 .Histo1D<float>("x");
1011
1012// request the generation of varied results from the nominal_hx
1013ROOT::RDF::Experimental::RResultMap<TH1D> hx = ROOT::RDF::Experimental::VariationsFor(nominal_hx);
1014
1015// the event loop runs here, upon first access to any of the results or varied results:
1016hx["nominal"].Draw(); // same effect as nominal_hx->Draw()
1017hx["pt:down"].Draw("SAME");
1018hx["pt:up"].Draw("SAME");
1019~~~
1020
1021A list of variation "tags" is passed as the last argument to Vary(). The tags give names to the varied values that are returned
1022as elements of an RVec of the appropriate C++ type. The number of variation tags must correspond to the number of elements of
1023this RVec (2 in the example above: the first element will correspond to the tag "down", the second
1024to the tag "up"). The _full_ variation name will be composed of the varied column name and the variation tags (e.g.
1025"pt:down", "pt:up" in this example). Python usage looks similar.
1026
1027Note how we use the "pt" column as usual in the Filter() and Define() calls and we simply use "x" as the value to fill
1028the resulting histogram. To produce the varied results, RDataFrame will automatically execute the Filter and Define
1029calls for each variation and fill the histogram with values and cuts that depend on the variation.
1030
1031There is no limitation to the complexity of a Vary() expression. Just like for the Define() and Filter() calls, users are
1032not limited to string expressions but they can also pass any valid C++ callable, including lambda functions and
1033complex functors. The callable can be applied to zero or more existing columns and it will always receive their
1034_nominal_ value in input.
1035
1036#### Varying multiple columns in lockstep
1037
1038In the following Python snippet we use the Vary() signature that allows varying multiple columns simultaneously or
1039"in lockstep":
1040
1041~~~{.python}
1042df.Vary(["pt", "eta"],
1043 "RVec<RVecF>{{pt*0.9, pt*1.1}, {eta*0.9, eta*1.1}}",
1044 variationTags=["down", "up"],
1045 variationName="ptAndEta")
1046~~~
1047
1048The expression returns an RVec of two RVecs: each inner vector contains the varied values for one column. The
1049inner vectors follow the same ordering as the column names that are passed as the first argument. Besides the variation tags, in
1050this case we also have to explicitly pass the variation name (here: "ptAndEta") as the default column name does not exist.
1051
1052The above call will produce variations "ptAndEta:down" and "ptAndEta:up".
1053
1054#### Combining multiple variations
1055
1056Even if a result depends on multiple variations, only one variation is applied at a time, i.e. there will be no result produced
1057by applying multiple systematic variations at the same time.
1058For example, in the following example snippet, the RResultMap instance `all_h` will contain keys "nominal", "pt:down",
1059"pt:up", "eta:0", "eta:1", but no "pt:up&&eta:0" or similar:
1060
1061~~~{.cpp}
1062auto df = _df.Vary("pt",
1063 "ROOT::RVecD{pt*0.9, pt*1.1}",
1064 {"down", "up"})
1065 .Vary("eta",
1066 [](float eta) { return RVecF{eta*0.9f, eta*1.1f}; },
1067 {"eta"},
1068 2);
1069
1070auto nom_h = df.Histo2D(histoModel, "pt", "eta");
1071auto all_hs = VariationsFor(nom_h);
1072all_hs.GetKeys(); // returns {"nominal", "pt:down", "pt:up", "eta:0", "eta:1"}
1073~~~
1074
1075Note how we passed the integer `2` instead of a list of variation tags to the second Vary() invocation: this is a
1076shorthand that automatically generates tags 0 to N-1 (in this case 0 and 1).
1077
1078\note Currently, VariationsFor() and RResultMap are in the `ROOT::RDF::Experimental` namespace, to indicate that these
1079 interfaces might still evolve and improve based on user feedback. We expect that some aspects of the related
1080 programming model will be streamlined in future versions.
1081
1082\note Currently, the results of a Snapshot(), Report() or Display() call cannot be varied (i.e. it is not possible to
1083 call \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" on them. These limitations will be lifted in future releases.
1084
1085See the Vary() method for more information and [this tutorial](https://root.cern/doc/master/df106__HiggsToFourLeptons_8C.html)
1086for an example usage of Vary and \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" in the analysis.
1087
1088\anchor rnode
1089### RDataFrame objects as function arguments and return values
1090RDataFrame variables/nodes are relatively cheap to copy and it's possible to both pass them to (or move them into)
1091functions and to return them from functions. However, in general each dataframe node will have a different C++ type,
1092which includes all available compile-time information about what that node does. One way to cope with this complication
1093is to use template functions and/or C++14 auto return types:
1094~~~{.cpp}
1095template <typename RDF>
1096auto ApplySomeFilters(RDF df)
1097{
1098 return df.Filter("x > 0").Filter([](int y) { return y < 0; }, {"y"});
1099}
1100~~~
1101
1102A possibly simpler, C++11-compatible alternative is to take advantage of the fact that any dataframe node can be
1103converted (implicitly or via an explicit cast) to the common type ROOT::RDF::RNode:
1104~~~{.cpp}
1105// a function that conditionally adds a Range to an RDataFrame node.
1106RNode MaybeAddRange(RNode df, bool mustAddRange)
1107{
1108 return mustAddRange ? df.Range(1) : df;
1109}
1110// use as :
1111ROOT::RDataFrame df(10);
1112auto maybeRangedDF = MaybeAddRange(df, true);
1113~~~
1114
1115The conversion to ROOT::RDF::RNode is cheap, but it will introduce an extra virtual call during the RDataFrame event
1116loop (in most cases, the resulting performance impact should be negligible). Python users can perform the conversion with the helper function `ROOT.RDF.AsRNode`.
1117
1118\anchor RDFCollections
1119### Storing RDataFrame objects in collections
1120
1121ROOT::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>`:
1122
1123~~~{.cpp}
1124std::vector<ROOT::RDF::RNode> dfs;
1125dfs.emplace_back(ROOT::RDataFrame(10));
1126dfs.emplace_back(dfs[0].Define("x", "42.f"));
1127~~~
1128
1129\anchor callbacks
1130### Executing callbacks every N events
1131It's possible to schedule execution of arbitrary functions (callbacks) during the event loop.
1132Callbacks can be used e.g. to inspect partial results of the analysis while the event loop is running,
1133drawing a partially-filled histogram every time a certain number of new entries is processed, or
1134displaying a progress bar while the event loop runs.
1135
1136For example one can draw an up-to-date version of a result histogram every 100 entries like this:
1137~~~{.cpp}
1138auto h = df.Histo1D("x");
1139TCanvas c("c","x hist");
1140h.OnPartialResult(100, [&c](TH1D &h_) { c.cd(); h_.Draw(); c.Update(); });
1141// event loop runs here, this final `Draw` is executed after the event loop is finished
1142h->Draw();
1143~~~
1144
1145Callbacks are registered to a ROOT::RDF::RResultPtr and must be callables that takes a reference to the result type as argument
1146and return nothing. RDataFrame will invoke registered callbacks passing partial action results as arguments to them
1147(e.g. a histogram filled with a part of the selected events).
1148
1149Read more on ROOT::RDF::RResultPtr::OnPartialResult() and ROOT::RDF::RResultPtr::OnPartialResultSlot().
1150
1151\anchor default-branches
1152### Default column lists
1153When constructing an RDataFrame object, it is possible to specify a **default column list** for your analysis, in the
1154usual form of a list of strings representing branch/column names. The default column list will be used as a fallback
1155whenever a list specific to the transformation/action is not present. RDataFrame will take as many of these columns as
1156needed, ignoring trailing extra names if present.
1157~~~{.cpp}
1158// use "b1" and "b2" as default columns
1159RDataFrame d1("myTree", "file.root", {"b1","b2"});
1160auto h = d1.Filter([](int b1, int b2) { return b1 > b2; }) // will act on "b1" and "b2"
1161 .Histo1D(); // will act on "b1"
1162
1163// just one default column this time
1164RDataFrame d2("myTree", "file.root", {"b1"});
1165auto min = d2.Filter([](double b2) { return b2 > 0; }, {"b2"}) // we can still specify non-default column lists
1166 .Min(); // returns the minimum value of "b1" for the filtered entries
1167~~~
1168
1169\anchor helper-cols
1170### Special helper columns: rdfentry_ and rdfslot_
1171Every instance of RDataFrame is created with two special columns called `rdfentry_` and `rdfslot_`. The `rdfentry_`
1172column is of type `ULong64_t` and it holds the current entry number while `rdfslot_` is an `unsigned int`
1173holding the index of the current data processing slot.
1174For backwards compatibility reasons, the names `tdfentry_` and `tdfslot_` are also accepted.
1175These columns are ignored by operations such as [Cache](classROOT_1_1RDF_1_1RInterface.html#aaaa0a7bb8eb21315d8daa08c3e25f6c9)
1176or [Snapshot](classROOT_1_1RDF_1_1RInterface.html#a233b7723e498967f4340705d2c4db7f8).
1177
1178\warning Note that in multi-thread event loops the values of `rdfentry_` _do not_ correspond to what would be the entry numbers
1179of a TChain constructed over the same set of ROOT files, as the entries are processed in an unspecified order.
1180
1181\anchor jitting
1182### Just-in-time compilation: column type inference and explicit declaration of column types
1183C++ is a statically typed language: all types must be known at compile-time. This includes the types of the TTree
1184branches we want to work on. For filters, defined columns and some of the actions, **column types are deduced from the
1185signature** of the relevant filter function/temporary column expression/action function:
1186~~~{.cpp}
1187// here b1 is deduced to be `int` and b2 to be `double`
1188df.Filter([](int x, double y) { return x > 0 && y < 0.; }, {"b1", "b2"});
1189~~~
1190If we specify an incorrect type for one of the columns, an exception with an informative message will be thrown at
1191runtime, when the column value is actually read from the dataset: RDataFrame detects type mismatches. The same would
1192happen if we swapped the order of "b1" and "b2" in the column list passed to Filter().
1193
1194Certain actions, on the other hand, do not take a function as argument (e.g. Histo1D()), so we cannot deduce the type of
1195the column at compile-time. In this case **RDataFrame infers the type of the column** from the TTree itself. This
1196is why we never needed to specify the column types for all actions in the above snippets.
1197
1198When the column type is not a common one such as `int`, `double`, `char` or `float` it is nonetheless good practice to
1199specify it as a template parameter to the action itself, like this:
1200~~~{.cpp}
1201df.Histo1D("b1"); // OK, the type of "b1" is deduced at runtime
1202df.Min<MyNumber_t>("myObject"); // OK, "myObject" is deduced to be of type `MyNumber_t`
1203~~~
1204
1205Deducing types at runtime requires the just-in-time compilation of the relevant actions, which has a small runtime
1206overhead, so specifying the type of the columns as template parameters to the action is good practice when performance is a goal.
1207
1208When strings are passed as expressions to Filter() or Define(), fundamental types are passed as constants. This avoids certaincommon mistakes such as typing `x = 0` rather than `x == 0`:
1209
1210~~~{.cpp}
1211// this throws an error (note the typo)
1212df.Define("x", "0").Filter("x = 0");
1213~~~
1214
1215\anchor generic-actions
1216### User-defined custom actions
1217RDataFrame strives to offer a comprehensive set of standard actions that can be performed on each event. At the same
1218time, it allows users to inject their own action code to perform arbitrarily complex data reductions.
1219
1220#### Implementing custom actions with Book()
1221
1222Through the Book() method, users can implement a custom action and have access to the same features
1223that built-in RDataFrame actions have, e.g. hooks to events related to the start, end and execution of the
1224event loop, or the possibility to return a lazy RResultPtr to an arbitrary type of result:
1225
1226~~~{.cpp}
1227#include <ROOT/RDataFrame.hxx>
1228#include <memory>
1229
1230class MyCounter : public ROOT::Detail::RDF::RActionImpl<MyCounter> {
1231 std::shared_ptr<int> fFinalResult = std::make_shared<int>(0);
1232 std::vector<int> fPerThreadResults;
1233
1234public:
1235 // We use a public type alias to advertise the type of the result of this action
1236 using Result_t = int;
1237
1238 MyCounter(unsigned int nSlots) : fPerThreadResults(nSlots) {}
1239
1240 // Called before the event loop to retrieve the address of the result that will be filled/generated.
1241 std::shared_ptr<int> GetResultPtr() const { return fFinalResult; }
1242
1243 // Called at the beginning of the event loop.
1244 void Initialize() {}
1245
1246 // Called at the beginning of each processing task.
1247 void InitTask(TTreeReader *, int) {}
1248
1249 /// Called at every entry.
1250 void Exec(unsigned int slot)
1251 {
1252 fPerThreadResults[slot]++;
1253 }
1254
1255 // Called at the end of the event loop.
1256 void Finalize()
1257 {
1258 *fFinalResult = std::accumulate(fPerThreadResults.begin(), fPerThreadResults.end(), 0);
1259 }
1260
1261 // Called by RDataFrame to retrieve the name of this action.
1262 std::string GetActionName() const { return "MyCounter"; }
1263};
1264
1265int main() {
1266 ROOT::RDataFrame df(10);
1267 ROOT::RDF::RResultPtr<int> resultPtr = df.Book<>(MyCounter{df.GetNSlots()}, {});
1268 // The GetValue call triggers the event loop
1269 std::cout << "Number of processed entries: " << resultPtr.GetValue() << std::endl;
1270}
1271~~~
1272
1273See the Book() method for more information and [this tutorial](https://root.cern/doc/master/df018__customActions_8C.html)
1274for a more complete example.
1275
1276#### Injecting arbitrary code in the event loop with Foreach() and ForeachSlot()
1277
1278Foreach() takes a callable (lambda expression, free function, functor...) and a list of columns and
1279executes the callable on the values of those columns for each event that passes all upstream selections.
1280It can be used to perform actions that are not already available in the interface. For example, the following snippet
1281evaluates the root mean square of column "x":
1282~~~{.cpp}
1283// Single-thread evaluation of RMS of column "x" using Foreach
1284double sumSq = 0.;
1285unsigned int n = 0;
1286df.Foreach([&sumSq, &n](double x) { ++n; sumSq += x*x; }, {"x"});
1287std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1288~~~
1289In multi-thread runs, users are responsible for the thread-safety of the expression passed to Foreach():
1290thread will execute the expression concurrently.
1291The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of `rms`; but
1292this is probably too much head-scratch for such a simple operation.
1293
1294ForeachSlot() can help in this situation. It is an alternative version of Foreach() for which the function takes an
1295additional "processing slot" parameter besides the columns it should be applied to. RDataFrame
1296guarantees that ForeachSlot() will invoke the user expression with different `slot` parameters for different concurrent
1297executions (see [Special helper columns: rdfentry_ and rdfslot_](\ref helper-cols) for more information on the slot parameter).
1298We can take advantage of ForeachSlot() to evaluate a thread-safe root mean square of column "x":
1299~~~{.cpp}
1300// Thread-safe evaluation of RMS of column "x" using ForeachSlot
1301ROOT::EnableImplicitMT();
1302const unsigned int nSlots = df.GetNSlots();
1303std::vector<double> sumSqs(nSlots, 0.);
1304std::vector<unsigned int> ns(nSlots, 0);
1305
1306df.ForeachSlot([&sumSqs, &ns](unsigned int slot, double x) { sumSqs[slot] += x*x; ns[slot] += 1; }, {"x"});
1307double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
1308unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
1309std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1310~~~
1311Notice how we created one `double` variable for each processing slot and later merged their results via `std::accumulate`.
1312
1313
1314\anchor friends
1315### Dataset joins with friend trees
1316
1317Vertically concatenating multiple trees that have the same columns (creating a logical dataset with the same columns and
1318more rows) is trivial in RDataFrame: just pass the tree name and a list of file names to RDataFrame's constructor, or create a TChain
1319out of the desired trees and pass that to RDataFrame.
1320
1321Horizontal concatenations of trees or chains (creating a logical dataset with the same number of rows and the union of the
1322columns of multiple trees) leverages TTree's "friend" mechanism.
1323
1324Simple joins of trees that do not have the same number of rows are also possible with indexed friend trees (see below).
1325
1326To use friend trees in RDataFrame, set up trees with the appropriate relationships and then instantiate an RDataFrame
1327with the main tree:
1328
1329~~~{.cpp}
1330TTree main([...]);
1331TTree friend([...]);
1332main.AddFriend(&friend, "myFriend");
1333
1334RDataFrame df(main);
1335auto df2 = df.Filter("myFriend.MyCol == 42");
1336~~~
1337
1338The same applies for TChains. Columns coming from the friend trees can be referred to by their full name, like in the example above,
1339or the friend tree name can be omitted in case the column name is not ambiguous (e.g. "MyCol" could be used instead of
1340"myFriend.MyCol" in the example above if there is no column "MyCol" in the main tree).
1341
1342\note A common source of confusion is that trees that are written out from a multi-thread Snapshot() call will have their
1343 entries (block-wise) shuffled with respect to the original tree. Such trees cannot be used as friends of the original
1344 one: rows will be mismatched.
1345
1346Indexed friend trees provide a way to perform simple joins of multiple trees over a common column.
1347When a certain entry in the main tree (or chain) is loaded, the friend trees (or chains) will then load an entry where the
1348"index" columns have a value identical to the one in the main one. For example, in Python:
1349
1350~~~{.py}
1351main_tree = ...
1352aux_tree = ...
1353
1354# If a friend tree has an index on `commonColumn`, when the main tree loads
1355# a given row, it also loads the row of the friend tree that has the same
1356# value of `commonColumn`
1357aux_tree.BuildIndex("commonColumn")
1358
1359mainTree.AddFriend(aux_tree)
1360
1361df = ROOT.RDataFrame(mainTree)
1362~~~
1363
1364RDataFrame supports indexed friend TTrees from ROOT v6.24 in single-thread mode and from v6.28/02 in multi-thread mode.
1365
1366\anchor other-file-formats
1367### Reading data formats other than ROOT trees
1368RDataFrame can be interfaced with RDataSources. The ROOT::RDF::RDataSource interface defines an API that RDataFrame can use to read arbitrary columnar data formats.
1369
1370RDataFrame calls into concrete RDataSource implementations to retrieve information about the data, retrieve (thread-local) readers or "cursors" for selected columns
1371and to advance the readers to the desired data entry.
1372Some predefined RDataSources are natively provided by ROOT such as the ROOT::RDF::RCsvDS which allows to read comma separated files:
1373~~~{.cpp}
1374auto tdf = ROOT::RDF::FromCSV("MuRun2010B.csv");
1375auto filteredEvents =
1376 tdf.Filter("Q1 * Q2 == -1")
1377 .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
1378auto h = filteredEvents.Histo1D("m");
1379h->Draw();
1380~~~
1381
1382See also FromNumpy (Python-only), FromRNTuple(), FromArrow(), FromSqlite().
1383
1384\anchor callgraphs
1385### Computation graphs (storing and reusing sets of transformations)
1386
1387As we saw, transformed dataframes can be stored as variables and reused multiple times to create modified versions of the dataset. This implicitly defines a **computation graph** in which
1388several paths of filtering/creation of columns are executed simultaneously, and finally aggregated results are produced.
1389
1390RDataFrame detects when several actions use the same filter or the same defined column, and **only evaluates each
1391filter or defined column once per event**, regardless of how many times that result is used down the computation graph.
1392Objects read from each column are **built once and never copied**, for maximum efficiency.
1393When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated,
1394so it might be advisable to put the strictest filters first in the graph.
1395
1396\anchor representgraph
1397### Visualizing the computation graph
1398It is possible to print the computation graph from any node to obtain a [DOT (graphviz)](https://en.wikipedia.org/wiki/DOT_(graph_description_language)) representation either on the standard output
1399or in a file.
1400
1401Invoking the function ROOT::RDF::SaveGraph() on any node that is not the head node, the computation graph of the branch
1402the node belongs to is printed. By using the head node, the entire computation graph is printed.
1403
1404Following there is an example of usage:
1405~~~{.cpp}
1406// First, a sample computational graph is built
1407ROOT::RDataFrame df("tree", "f.root");
1408
1409auto df2 = df.Define("x", []() { return 1; })
1410 .Filter("col0 % 1 == col0")
1411 .Filter([](int b1) { return b1 <2; }, {"cut1"})
1412 .Define("y", []() { return 1; });
1413
1414auto count = df2.Count();
1415
1416// Prints the graph to the rd1.dot file in the current directory
1417ROOT::RDF::SaveGraph(df, "./mydot.dot");
1418// Prints the graph to standard output
1419ROOT::RDF::SaveGraph(df);
1420~~~
1421
1422The generated graph can be rendered using one of the graphviz filters, e.g. `dot`. For instance, the image below can be generated with the following command:
1423~~~{.sh}
1424$ dot -Tpng computation_graph.dot -ocomputation_graph.png
1425~~~
1426
1427\image html RDF_Graph2.png
1428
1429\anchor rdf-logging
1430### Activating RDataFrame execution logs
1431
1432RDataFrame has experimental support for verbose logging of the event loop runtimes and other interesting related information. It is activated as follows:
1433~~~{.cpp}
1434#include <ROOT/RLogger.hxx>
1435
1436// this increases RDF's verbosity level as long as the `verbosity` variable is in scope
1437auto verbosity = ROOT::Experimental::RLogScopedVerbosity(ROOT::Detail::RDF::RDFLogChannel(), ROOT::Experimental::ELogLevel::kInfo);
1438~~~
1439
1440or in Python:
1441~~~{.python}
1442import ROOT
1443
1444verbosity = ROOT.Experimental.RLogScopedVerbosity(ROOT.Detail.RDF.RDFLogChannel(), ROOT.Experimental.ELogLevel.kInfo)
1445~~~
1446
1447More information (e.g. start and end of each multi-thread task) is printed using `ELogLevel.kDebug` and even more
1448(e.g. a full dump of the generated code that RDataFrame just-in-time-compiles) using `ELogLevel.kDebug+10`.
1449
1450\anchor rdf-from-spec
1451### Creating an RDataFrame from a dataset specification file
1452
1453RDataFrame can be created using a dataset specification JSON file:
1454
1455~~~{.python}
1456import ROOT
1457
1458df = ROOT.RDF.Experimental.FromSpec("spec.json")
1459~~~
1460
1461The input dataset specification JSON file needs to be provided by the user and it describes all necessary samples and
1462their associated metadata information. The main required key is the "samples" (at least one sample is needed) and the
1463required sub-keys for each sample are "trees" and "files". Additionally, one can specify a metadata dictionary for each
1464sample in the "metadata" key.
1465
1466A simple example for the formatting of the specification in the JSON file is the following:
1467
1468~~~{.cpp}
1469{
1470 "samples": {
1471 "sampleA": {
1472 "trees": ["tree1", "tree2"],
1473 "files": ["file1.root", "file2.root"],
1474 "metadata": {
1475 "lumi": 10000.0,
1476 "xsec": 1.0,
1477 "sample_category" = "data"
1478 }
1479 },
1480 "sampleB": {
1481 "trees": ["tree3", "tree4"],
1482 "files": ["file3.root", "file4.root"],
1483 "metadata": {
1484 "lumi": 0.5,
1485 "xsec": 1.5,
1486 "sample_category" = "MC_background"
1487 }
1488 }
1489 }
1490}
1491~~~
1492
1493The metadata information from the specification file can be then accessed using the DefinePerSample function.
1494For example, to access luminosity information (stored as a double):
1495
1496~~~{.python}
1497df.DefinePerSample("lumi", 'rdfsampleinfo_.GetD("lumi")')
1498~~~
1499
1500or sample_category information (stored as a string):
1501
1502~~~{.python}
1503df.DefinePerSample("sample_category", 'rdfsampleinfo_.GetS("sample_category")')
1504~~~
1505
1506or directly the filename:
1507
1508~~~{.python}
1509df.DefinePerSample("name", "rdfsampleinfo_.GetSampleName()")
1510~~~
1511
1512An example implementation of the "FromSpec" method is available in tutorial: df106_HiggstoFourLeptons.py, which also
1513provides a corresponding exemplary JSON file for the dataset specification.
1514
1515\anchor progressbar
1516### Adding a progress bar
1517
1518A progress bar showing the processed event statistics can be added to any RDataFrame program.
1519The event statistics include elapsed time, currently processed file, currently processed events, the rate of event processing
1520and an estimated remaining time (per file being processed). It is recorded and printed in the terminal every m events and every
1521n seconds (by default m = 1000 and n = 1). The ProgressBar can be also added when the multithread (MT) mode is enabled.
1522
1523ProgressBar is added after creating the dataframe object (df):
1524~~~{.cpp}
1525ROOT::RDataFrame df("tree", "file.root");
1526ROOT::RDF::Experimental::AddProgressBar(df);
1527~~~
1528
1529Alternatively, RDataFrame can be cast to an RNode first, giving the user more flexibility
1530For example, it can be called at any computational node, such as Filter or Define, not only the head node,
1531with no change to the ProgressBar function itself (please see the [Efficient analysis in Python](#python)
1532section for appropriate usage in Python):
1533~~~{.cpp}
1534ROOT::RDataFrame df("tree", "file.root");
1535auto df_1 = ROOT::RDF::RNode(df.Filter("x>1"));
1536ROOT::RDF::Experimental::AddProgressBar(df_1);
1537~~~
1538Examples of implemented progress bars can be seen by running [Higgs to Four Lepton tutorial](https://root.cern/doc/master/df106__HiggsToFourLeptons_8py_source.html) and [Dimuon tutorial](https://root.cern/doc/master/df102__NanoAODDimuonAnalysis_8C.html).
1539
1540*/
1541// clang-format on
1542
1543namespace ROOT {
1544
1546using ColumnNamesPtr_t = std::shared_ptr<const ColumnNames_t>;
1547
1548////////////////////////////////////////////////////////////////////////////
1549/// \brief Build the dataframe.
1550/// \param[in] treeName Name of the tree contained in the directory
1551/// \param[in] dirPtr TDirectory where the tree is stored, e.g. a TFile.
1552/// \param[in] defaultColumns Collection of default columns.
1553///
1554/// The default columns are looked at in case no column is specified in the
1555/// booking of actions or transformations.
1556/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1557RDataFrame::RDataFrame(std::string_view treeName, TDirectory *dirPtr, const ColumnNames_t &defaultColumns)
1558 : RInterface(std::make_shared<RDFDetail::RLoopManager>(nullptr, defaultColumns))
1559{
1560 if (!dirPtr) {
1561 auto msg = "Invalid TDirectory!";
1562 throw std::runtime_error(msg);
1563 }
1564 const std::string treeNameInt(treeName);
1565 auto tree = static_cast<TTree *>(dirPtr->Get(treeNameInt.c_str()));
1566 if (!tree) {
1567 auto msg = "Tree \"" + treeNameInt + "\" cannot be found!";
1568 throw std::runtime_error(msg);
1569 }
1570 GetProxiedPtr()->SetTree(std::shared_ptr<TTree>(tree, [](TTree *) {}));
1571}
1572
1573////////////////////////////////////////////////////////////////////////////
1574/// \brief Build the dataframe.
1575/// \param[in] treeName Name of the tree contained in the directory
1576/// \param[in] filenameglob TDirectory where the tree is stored, e.g. a TFile.
1577/// \param[in] defaultColumns Collection of default columns.
1578///
1579/// The filename glob supports the same type of expressions as TChain::Add(), and it is passed as-is to TChain's
1580/// constructor.
1581///
1582/// The default columns are looked at in case no column is specified in the
1583/// booking of actions or transformations.
1584/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1585#ifdef R__HAS_ROOT7
1586RDataFrame::RDataFrame(std::string_view treeName, std::string_view fileNameGlob, const ColumnNames_t &defaultColumns)
1587 : RInterface(ROOT::Detail::RDF::CreateLMFromFile(treeName, fileNameGlob, defaultColumns))
1588{
1589}
1590#else
1591RDataFrame::RDataFrame(std::string_view treeName, std::string_view fileNameGlob, const ColumnNames_t &defaultColumns)
1592 : RInterface(ROOT::Detail::RDF::CreateLMFromTTree(treeName, fileNameGlob, defaultColumns))
1593{
1594}
1595#endif
1596
1597////////////////////////////////////////////////////////////////////////////
1598/// \brief Build the dataframe.
1599/// \param[in] treeName Name of the tree contained in the directory
1600/// \param[in] fileglobs Collection of file names of filename globs
1601/// \param[in] defaultColumns Collection of default columns.
1602///
1603/// The filename globs support the same type of expressions as TChain::Add(), and each glob is passed as-is
1604/// to TChain's constructor.
1605///
1606/// The default columns are looked at in case no column is specified in the booking of actions or transformations.
1607/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1608#ifdef R__HAS_ROOT7
1609RDataFrame::RDataFrame(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1610 const ColumnNames_t &defaultColumns)
1611 : RInterface(ROOT::Detail::RDF::CreateLMFromFile(datasetName, fileNameGlobs, defaultColumns))
1612{
1613}
1614#else
1615RDataFrame::RDataFrame(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1616 const ColumnNames_t &defaultColumns)
1617 : RInterface(ROOT::Detail::RDF::CreateLMFromTTree(datasetName, fileNameGlobs, defaultColumns))
1618{
1619}
1620#endif
1621
1622////////////////////////////////////////////////////////////////////////////
1623/// \brief Build the dataframe.
1624/// \param[in] tree The tree or chain to be studied.
1625/// \param[in] defaultColumns Collection of default column names to fall back to when none is specified.
1626///
1627/// The default columns are looked at in case no column is specified in the
1628/// booking of actions or transformations.
1629/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1630RDataFrame::RDataFrame(TTree &tree, const ColumnNames_t &defaultColumns)
1631 : RInterface(std::make_shared<RDFDetail::RLoopManager>(&tree, defaultColumns))
1632{
1633}
1634
1635//////////////////////////////////////////////////////////////////////////
1636/// \brief Build a dataframe that generates numEntries entries.
1637/// \param[in] numEntries The number of entries to generate.
1638///
1639/// An empty-source dataframe constructed with a number of entries will
1640/// generate those entries on the fly when some action is triggered,
1641/// and it will do so for all the previously-defined columns.
1642/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1644 : RInterface(std::make_shared<RDFDetail::RLoopManager>(numEntries))
1645
1646{
1647}
1648
1649//////////////////////////////////////////////////////////////////////////
1650/// \brief Build dataframe associated to data source.
1651/// \param[in] ds The data source object.
1652/// \param[in] defaultColumns Collection of default column names to fall back to when none is specified.
1653///
1654/// A dataframe associated to a data source will query it to access column values.
1655/// \see ROOT::RDF::RInterface for the documentation of the methods available.
1656RDataFrame::RDataFrame(std::unique_ptr<ROOT::RDF::RDataSource> ds, const ColumnNames_t &defaultColumns)
1657 : RInterface(std::make_shared<RDFDetail::RLoopManager>(std::move(ds), defaultColumns))
1658{
1659}
1660
1661//////////////////////////////////////////////////////////////////////////
1662/// \brief Build dataframe from an RDatasetSpec object.
1663/// \param[in] spec The dataset specification object.
1664///
1665/// A dataset specification includes trees and file names,
1666/// as well as an optional friend list and/or entry range.
1667///
1668/// ### Example usage from Python:
1669/// ~~~{.py}
1670/// spec = (
1671/// ROOT.RDF.Experimental.RDatasetSpec()
1672/// .AddSample(("data", "tree", "file.root"))
1673/// .WithGlobalFriends("friendTree", "friend.root", "alias")
1674/// .WithGlobalRange((100, 200))
1675/// )
1676/// df = ROOT.RDataFrame(spec)
1677/// ~~~
1678///
1679/// See also ROOT::RDataFrame::FromSpec().
1681 : RInterface(std::make_shared<RDFDetail::RLoopManager>(std::move(spec)))
1682{
1683}
1684
1686{
1687 // If any node of the computation graph associated with this RDataFrame
1688 // declared code to jit, we need to make sure the compilation actually
1689 // happens. For example, a jitted Define could have been booked but
1690 // if the computation graph is not actually run then the code of the
1691 // Define node is not jitted. This in turn would cause memory leaks.
1692 // See https://github.com/root-project/root/issues/15399
1693 fLoopManager->Jit();
1694}
1695
1696namespace RDF {
1697namespace Experimental {
1698
1699////////////////////////////////////////////////////////////////////////////
1700/// \brief Create the RDataFrame from the dataset specification file.
1701/// \param[in] jsonFile Path to the dataset specification JSON file.
1702///
1703/// The input dataset specification JSON file must include a number of keys that
1704/// describe all the necessary samples and their associated metadata information.
1705/// The main key, "samples", is required and at least one sample is needed. Each
1706/// sample must have at least one key "trees" and at least one key "files" from
1707/// which the data is read. Optionally, one or more metadata information can be
1708/// added, as well as the friend list information.
1709///
1710/// ### Example specification file JSON:
1711/// The following is an example of the dataset specification JSON file formatting:
1712///~~~{.cpp}
1713/// {
1714/// "samples": {
1715/// "sampleA": {
1716/// "trees": ["tree1", "tree2"],
1717/// "files": ["file1.root", "file2.root"],
1718/// "metadata": {"lumi": 1.0, }
1719/// },
1720/// "sampleB": {
1721/// "trees": ["tree3", "tree4"],
1722/// "files": ["file3.root", "file4.root"],
1723/// "metadata": {"lumi": 0.5, }
1724/// },
1725/// ...
1726/// },
1727/// }
1728///~~~
1729ROOT::RDataFrame FromSpec(const std::string &jsonFile)
1730{
1731 const nlohmann::ordered_json fullData = nlohmann::ordered_json::parse(std::ifstream(jsonFile));
1732 if (!fullData.contains("samples") || fullData["samples"].empty()) {
1733 throw std::runtime_error(
1734 R"(The input specification does not contain any samples. Please provide the samples in the specification like:
1735{
1736 "samples": {
1737 "sampleA": {
1738 "trees": ["tree1", "tree2"],
1739 "files": ["file1.root", "file2.root"],
1740 "metadata": {"lumi": 1.0, }
1741 },
1742 "sampleB": {
1743 "trees": ["tree3", "tree4"],
1744 "files": ["file3.root", "file4.root"],
1745 "metadata": {"lumi": 0.5, }
1746 },
1747 ...
1748 },
1749})");
1750 }
1751
1752 RDatasetSpec spec;
1753 for (const auto &keyValue : fullData["samples"].items()) {
1754 const std::string &sampleName = keyValue.key();
1755 const auto &sample = keyValue.value();
1756 // TODO: if requested in https://github.com/root-project/root/issues/11624
1757 // allow union-like types for trees and files, see: https://github.com/nlohmann/json/discussions/3815
1758 if (!sample.contains("trees")) {
1759 throw std::runtime_error("A list of tree names must be provided for sample " + sampleName + ".");
1760 }
1761 std::vector<std::string> trees = sample["trees"];
1762 if (!sample.contains("files")) {
1763 throw std::runtime_error("A list of files must be provided for sample " + sampleName + ".");
1764 }
1765 std::vector<std::string> files = sample["files"];
1766 if (!sample.contains("metadata")) {
1767 spec.AddSample(RSample{sampleName, trees, files});
1768 } else {
1769 RMetaData m;
1770 for (const auto &metadata : sample["metadata"].items()) {
1771 const auto &val = metadata.value();
1772 if (val.is_string())
1773 m.Add(metadata.key(), val.get<std::string>());
1774 else if (val.is_number_integer())
1775 m.Add(metadata.key(), val.get<int>());
1776 else if (val.is_number_float())
1777 m.Add(metadata.key(), val.get<double>());
1778 else
1779 throw std::logic_error("The metadata keys can only be of type [string|int|double].");
1780 }
1781 spec.AddSample(RSample{sampleName, trees, files, m});
1782 }
1783 }
1784 if (fullData.contains("friends")) {
1785 for (const auto &friends : fullData["friends"].items()) {
1786 std::string alias = friends.key();
1787 std::vector<std::string> trees = friends.value()["trees"];
1788 std::vector<std::string> files = friends.value()["files"];
1789 if (files.size() != trees.size() && trees.size() > 1)
1790 throw std::runtime_error("Mismatch between trees and files in a friend.");
1791 spec.WithGlobalFriends(trees, files, alias);
1792 }
1793 }
1794
1795 if (fullData.contains("range")) {
1796 std::vector<int> range = fullData["range"];
1797
1798 if (range.size() == 1)
1799 spec.WithGlobalRange({range[0]});
1800 else if (range.size() == 2)
1801 spec.WithGlobalRange({range[0], range[1]});
1802 }
1803 return ROOT::RDataFrame(spec);
1804}
1805
1806} // namespace Experimental
1807} // namespace RDF
1808
1809} // namespace ROOT
1810
1811namespace cling {
1812//////////////////////////////////////////////////////////////////////////
1813/// Print an RDataFrame at the prompt
1814std::string printValue(ROOT::RDataFrame *df)
1815{
1816 // The loop manager is never null, except when its construction failed.
1817 // This can happen e.g. if the constructor of RLoopManager that expects
1818 // a file name is used and that file doesn't exist. This point is usually
1819 // not even reached in that situation, since the exception thrown by the
1820 // constructor will also stop execution of the program. But it can still
1821 // be reached at the prompt, if the user tries to print the RDataFrame
1822 // variable after an incomplete initialization.
1823 auto *lm = df->GetLoopManager();
1824 if (!lm) {
1825 throw std::runtime_error("Cannot print information about this RDataFrame, "
1826 "it was not properly created. It must be discarded.");
1827 }
1828 auto *tree = lm->GetTree();
1829 auto defCols = lm->GetDefaultColumnNames();
1830
1831 std::ostringstream ret;
1832 if (tree) {
1833 ret << "A data frame built on top of the " << tree->GetName() << " dataset.";
1834 if (!defCols.empty()) {
1835 if (defCols.size() == 1)
1836 ret << "\nDefault column: " << defCols[0];
1837 else {
1838 ret << "\nDefault columns:\n";
1839 for (auto &&col : defCols) {
1840 ret << " - " << col << "\n";
1841 }
1842 }
1843 }
1844 } else if (auto ds = df->fDataSource) {
1845 ret << "A data frame associated to the data source \"" << cling::printValue(ds) << "\"";
1846 } else {
1847 ret << "An empty data frame that will create " << lm->GetNEmptyEntries() << " entries\n";
1848 }
1849
1850 return ret.str();
1851}
1852} // namespace cling
unsigned long long ULong64_t
Definition RtypesCore.h:81
The head node of a RDF computation graph.
The dataset specification for RDataFrame.
RDatasetSpec & WithGlobalFriends(const std::string &treeName, const std::string &fileNameGlob, const std::string &alias="")
Add friend tree to RDatasetSpec object.
RDatasetSpec & AddSample(RSample sample)
Add sample (RSample class object) to the RDatasetSpec object.
RDatasetSpec & WithGlobalRange(const RDatasetSpec::REntryRange &entryRange={})
Create an RDatasetSpec object for a given range of entries.
Class behaving as a heterogenuous dictionary to store the metadata of a dataset.
Definition RMetaData.hxx:50
Class representing a sample which is a grouping of trees and their fileglobs, and,...
Definition RSample.hxx:39
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > fLoopManager
< The RLoopManager at the root of this computation graph. Never null.
RDataSource * fDataSource
Non-owning pointer to a data-source object. Null if no data-source. RLoopManager has ownership of the...
RDFDetail::RLoopManager * GetLoopManager() const
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
RDataFrame(std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultColumns={})
Build the dataframe.
ROOT::RDF::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
ROOT::RDataFrame FromSpec(const std::string &jsonFile)
Factory method to create an RDataFrame from a JSON specification file.
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
TMarker m
Definition textangle.C:8