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