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