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