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