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