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
536ROOT provides convenience utility functions to work with RVec which are available in the [ROOT::VecOps namespace](\ref vecops)
537
538\anchor transformations
539## Transformations: manipulating data
540\anchor Filters
541### Filters
542A filter is created through a call to `Filter(f, columnList)` or `Filter(filterString)`. In the first overload, `f` can
543be a function, a lambda expression, a functor class, or any other callable object. It must return a `bool` signalling
544whether the event has passed the selection (`true`) or not (`false`). It should perform "read-only" operations on the
545columns, and should not have side-effects (e.g. modification of an external or static variable) to ensure correctness
546when implicit multi-threading is active. The second overload takes a string with a valid C++ expression in which column
547names are used as variable names (e.g. `Filter("x[0] + x[1] > 0")`). This is a convenience feature that comes with a
548certain runtime overhead: C++ code has to be generated on the fly from this expression before using it in the event
549loop. See the paragraph about "Just-in-time compilation" below for more information.
550
551RDataFrame only evaluates filters when necessary: if multiple filters are chained one after another, they are executed
552in order and the first one returning `false` causes the event to be discarded and triggers the processing of the next
553entry. If multiple actions or transformations depend on the same filter, that filter is not executed multiple times for
554each entry: after the first access it simply serves a cached result.
555
556\anchor named-filters-and-cutflow-reports
557#### Named filters and cutflow reports
558An optional string parameter `name` can be passed to the Filter() method to create a **named filter**. Named filters
559work as usual, but also keep track of how many entries they accept and reject.
560
561Statistics are retrieved through a call to the Report() method:
562
563- when Report() is called on the main RDataFrame object, it returns a ROOT::RDF::RResultPtr<RCutFlowReport> relative to all
564named filters declared up to that point
565- when called on a specific node (e.g. the result of a Define() or Filter()), it returns a ROOT::RDF::RResultPtr<RCutFlowReport>
566relative all named filters in the section of the chain between the main RDataFrame and that node (included).
567
568Stats are stored in the same order as named filters have been added to the graph, and *refer to the latest event-loop*
569that has been run using the relevant RDataFrame.
570
571\anchor ranges
572### Ranges
573When RDataFrame is not being used in a multi-thread environment (i.e. no call to EnableImplicitMT() was made),
574Range() transformations are available. These act very much like filters but instead of basing their decision on
575a filter expression, they rely on `begin`,`end` and `stride` parameters.
576
577- `begin`: initial entry number considered for this range.
578- `end`: final entry number (excluded) considered for this range. 0 means that the range goes until the end of the dataset.
579- `stride`: process one entry of the [begin, end) range every `stride` entries. Must be strictly greater than 0.
580
581The actual number of entries processed downstream of a Range() node will be `(end - begin)/stride` (or less if less
582entries than that are available).
583
584Note that ranges act "locally", not based on the global entry count: `Range(10,50)` means "skip the first 10 entries
585*that reach this node*, let the next 40 entries pass, then stop processing". If a range node hangs from a filter node,
586and the range has a `begin` parameter of 10, that means the range will skip the first 10 entries *that pass the
587preceding filter*.
588
589Ranges allow "early quitting": if all branches of execution of a functional graph reached their `end` value of
590processed entries, the event-loop is immediately interrupted. This is useful for debugging and quick data explorations.
591
592\anchor custom-columns
593### Custom columns
594Custom columns are created by invoking `Define(name, f, columnList)`. As usual, `f` can be any callable object
595(function, lambda expression, functor class...); it takes the values of the columns listed in `columnList` (a list of
596strings) as parameters, in the same order as they are listed in `columnList`. `f` must return the value that will be
597assigned to the temporary column.
598
599A new variable is created called `name`, accessible as if it was contained in the dataset from subsequent
600transformations/actions.
601
602Use cases include:
603- caching the results of complex calculations for easy and efficient multiple access
604- extraction of quantities of interest from complex objects
605- branch aliasing, i.e. changing the name of a branch
606
607An exception is thrown if the `name` of the new column/branch is already in use for another branch in the TTree.
608
609It is also possible to specify the quantity to be stored in the new temporary column as a C++ expression with the method
610`Define(name, expression)`. For example this invocation
611
612~~~{.cpp}
613df.Define("pt", "sqrt(px*px + py*py)");
614~~~
615
616will create a new column called "pt" the value of which is calculated starting from the columns px and py. The system
617builds a just-in-time compiled function starting from the expression after having deduced the list of necessary branches
618from the names of the variables specified by the user.
619
620#### Custom columns as function of slot and entry number
621
622It is possible to create custom columns also as a function of the processing slot and entry numbers. The methods that can
623be invoked are:
624- `DefineSlot(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, T1, T2, ...)`: the
625first parameter is the slot number which ranges from 0 to ROOT::GetThreadPoolSize() - 1.
626- `DefineSlotEntry(name, f, columnList)`. In this case the callable f has this signature `R(unsigned int, ULong64_t,
627T1, T2, ...)`: the first parameter is the slot number while the second one the number of the entry being processed.
628
629\anchor actions
630## Actions: getting results
631### Instant and lazy actions
632Actions can be **instant** or **lazy**. Instant actions are executed as soon as they are called, while lazy actions are
633executed whenever the object they return is accessed for the first time. As a rule of thumb, actions with a return value
634are lazy, the others are instant.
635
636### Return type of a lazy action
637
638When a lazy action is called, it returns a \link ROOT::RDF::RResultPtr ROOT::RDF::RResultPtr<T>\endlink, where T is the
639type of the result of the action. The final result will be stored in the `RResultPtr`, and can be retrieved by
640dereferencing it or via its `GetValue` method. Retrieving the result also starts the event loop if the result
641hasn't been produced yet.
642
643The RResultPtr shares ownership of the result object. To directly access result, use:
644~~~{.cpp}
645ROOT::RDF::RResultPtr<TH1D> histo = rdf.Histo1D(...);
646histo->Draw(); // Starts running the event loop
647~~~
648
649To return results from functions, a copy of the underlying shared_ptr can be obtained:
650~~~{.cpp}
651std::shared_ptr<TH1D> ProduceResult(const char *columnname) {
652 ROOT::RDF::RResultPtr<TH1D> histo = rdf.Histo1D(*h, columname);
653 return histo.GetSharedPtr(); // Runs the event loop
654}
655~~~
656If the result had been returned by reference or bare pointer, it would have gotten destroyed
657when the function exits.
658
659To share ownership but not produce the result ("keep it lazy"), copy the RResultPtr:
660~~~{.cpp}
661std::vector<RResultPtr<TH1D>> allHistograms;
662ROOT::RDF::RResultPtr<TH1D> BookHistogram(const char *columnname) {
663 ROOT::RDF::RResultPtr<TH1D> histo = rdf.Histo1D(*h, columname);
664 allHistograms.push_back(histo); // Will not produce the result yet
665 return histo;
666}
667~~~
668
669
670### Actions that return collections
671
672If the type of the return value of an action is a collection, e.g. `std::vector<int>`, you can iterate its elements
673directly through the wrapping `RResultPtr`:
674
675~~~{.cpp}
676ROOT::RDataFrame df{5};
677auto df1 = df.Define("x", []{ return 42; });
678for (const auto &el: df1.Take<int>("x")){
679 std::cout << "Element: " << el << "\n";
680}
681~~~
682
683~~~{.py}
684df = ROOT.RDataFrame(5).Define("x", "42")
685for el in df.Take[int]("x"):
686 print(f"Element: {el}")
687~~~
688
689### Actions and readers
690
691An action that needs values for its computations will request it from a reader, e.g. a column created via `Define` or
692available from the input dataset. The action will request values from each column of the list of input columns (either
693inferred or specified by the user), in order. For example:
694
695~~~{.cpp}
696ROOT::RDataFrame df{1};
697auto df1 = df.Define("x", []{ return 11; });
698auto df2 = df1.Define("y", []{ return 22; });
699auto graph = df2.Graph<int, int>("x","y");
700~~~
701
702The `Graph` action is going to request first the value from column "x", then that of column "y". Specifically, the order
703of execution of the operations of nodes in this branch of the computation graph is guaranteed to be top to bottom.
704
705\anchor rdf_distrdf
706## Distributed execution
707
708RDataFrame applications can be executed in parallel through distributed computing frameworks on a set of remote machines
709thanks to the Python package `ROOT.RDF.Distributed`. This **Python-only** package allows to scale the
710optimized performance RDataFrame can achieve on a single machine to multiple nodes at the same time. It is designed so
711that different backends can be easily plugged in, currently supporting [Apache Spark](http://spark.apache.org/) and
712[Dask](https://dask.org/). Here is a minimal example usage of distributed RDataFrame:
713
714~~~{.py}
715import ROOT
716from distributed import Client
717# It still accepts the same constructor arguments as traditional RDataFrame
718# but needs a client object which allows connecting to one of the supported
719# schedulers (read more info below)
720client = Client(...)
721df = ROOT.RDataFrame("mytree", "myfile.root", executor=client)
722
723# Continue the application with the traditional RDataFrame API
724sum = df.Filter("x > 10").Sum("y")
725h = df.Histo1D(("name", "title", 10, 0, 10), "x")
726
727print(sum.GetValue())
728h.Draw()
729~~~
730
731The main goal of this package is to support running any RDataFrame application distributedly. Nonetheless, not all
732parts of the RDataFrame API currently work with this package. The subset that is currently available is:
733- Alias
734- AsNumpy
735- Count
736- DefaultValueFor
737- Define
738- DefinePerSample
739- Filter
740- FilterAvailable
741- FilterMissing
742- Graph
743- Histo[1,2,3]D
744- HistoND, HistoNSparseD
745- Max
746- Mean
747- Min
748- Profile[1,2,3]D
749- Redefine
750- Snapshot
751- Stats
752- StdDev
753- Sum
754- Systematic variations: Vary and [VariationsFor](\ref ROOT::RDF::Experimental::VariationsFor).
755- Parallel submission of distributed graphs: [RunGraphs](\ref ROOT::RDF::RunGraphs).
756- Information about the dataframe: GetColumnNames.
757
758with support for more operations coming in the future. Currently, to the supported data sources belong TTree, TChain, RNTuple and RDatasetSpec.
759
760### Connecting to a Spark cluster
761
762In order to distribute the RDataFrame workload, you can connect to a Spark cluster you have access to through the
763official [Spark API](https://spark.apache.org/docs/latest/rdd-programming-guide.html#initializing-spark), then hook the
764connection instance to the distributed `RDataFrame` object like so:
765
766~~~{.py}
767import pyspark
768import ROOT
769
770# Create a SparkContext object with the right configuration for your Spark cluster
771conf = SparkConf().setAppName(appName).setMaster(master)
772sc = SparkContext(conf=conf)
773
774# The Spark RDataFrame constructor accepts an optional "sparkcontext" parameter
775# and it will distribute the application to the connected cluster
776df = ROOT.RDataFrame("mytree", "myfile.root", executor = sc)
777~~~
778
779Note that with the usage above the case of `executor = None` is not supported. One
780can explicitly create a `ROOT.RDF.Distributed.Spark.RDataFrame` object
781in order to get a default instance of
782[SparkContext](https://spark.apache.org/docs/latest/api/python/reference/api/pyspark.SparkContext.html)
783in case it is not already provided as argument.
784
785### Connecting to a Dask cluster
786
787Similarly, you can connect to a Dask cluster by creating your own connection object which internally operates with one
788of the cluster schedulers supported by Dask (more information in the
789[Dask distributed docs](http://distributed.dask.org/en/stable/)):
790
791~~~{.py}
792import ROOT
793from dask.distributed import Client
794# In a Python script the Dask client needs to be initalized in a context
795# Jupyter notebooks / Python session don't need this
796if __name__ == "__main__":
797 # With an already setup cluster that exposes a Dask scheduler endpoint
798 client = Client("dask_scheduler.domain.com:8786")
799
800 # The Dask RDataFrame constructor accepts the Dask Client object as an optional argument
801 df = ROOT.RDataFrame("mytree","myfile.root", executor=client)
802 # Proceed as usual
803 df.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
804~~~
805
806Note that with the usage above the case of `executor = None` is not supported. One
807can explicitly create a `ROOT.RDF.Distributed.Dask.RDataFrame` object
808in order to get a default instance of
809[distributed.Client](http://distributed.dask.org/en/stable/api.html#distributed.Client)
810in case it is not already provided as argument. This will run multiple processes
811on the local machine using all available cores.
812
813### Choosing the number of distributed tasks
814
815A distributed RDataFrame has internal logic to define in how many chunks the input dataset will be split before sending
816tasks to the distributed backend. Each task reads and processes one of said chunks. The logic is backend-dependent, but
817generically tries to infer how many cores are available in the cluster through the connection object. The number of
818tasks will be equal to the inferred number of cores. There are cases where the connection object of the chosen backend
819doesn't have information about the actual resources of the cluster. An example of this is when using Dask to connect to
820a batch system. The client object created at the beginning of the application does not automatically know how many cores
821will be available during distributed execution, since the jobs are submitted to the batch system after the creation of
822the connection. In such cases, the logic is to default to process the whole dataset in 2 tasks.
823
824The number of tasks submitted for distributed execution can be also set programmatically, by providing the optional
825keyword argument `npartitions` when creating the RDataFrame object. This parameter is accepted irrespectively of the
826backend used:
827
828~~~{.py}
829import ROOT
830
831if __name__ == "__main__":
832 # The `npartitions` optional argument tells the RDataFrame how many tasks are desired
833 df = ROOT.RDataFrame("mytree", "myfile.root", executor=SupportedExecutor(...), npartitions=NPARTITIONS)
834 # Proceed as usual
835 df.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
836~~~
837
838Note that when processing a TTree or TChain dataset, the `npartitions` value should not exceed the number of clusters in
839the dataset. The number of clusters in a TTree can be retrieved by typing `rootls -lt myfile.root` at a command line.
840
841### Distributed FromSpec
842
843RDataFrame 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
844jsonFile (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.
845If 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.
846For more information on FromSpec functionality itself please refer to [FromSpec](\ref rdf-from-spec) documentation. Note that adding metadata and friend information is supported,
847but adding the global range will not be respected in the distributed execution.
848
849Using spark:
850~~~{.py}
851import pyspark
852import ROOT
853
854conf = SparkConf().setAppName(appName).setMaster(master)
855sc = SparkContext(conf=conf)
856
857# The FromSpec function accepts an optional "sparkcontext" parameter
858# and it will distribute the application to the connected cluster
859df_fromspec = ROOT.RDF.Experimental.FromSpec("myspec.json", executor = sc)
860# Proceed as usual
861df_fromspec.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
862~~~
863
864Using dask:
865~~~{.py}
866import ROOT
867from dask.distributed import Client
868
869if __name__ == "__main__":
870 client = Client("dask_scheduler.domain.com:8786")
871
872 # The FromSpec function accepts the Dask Client object as an optional argument
873 df_fromspec = ROOT.RDF.Experimental.FromSpec("myspec.json", executor=client)
874 # Proceed as usual
875 df_fromspec.Define("x","someoperation").Histo1D(("name", "title", 10, 0, 10), "x")
876~~~
877
878### Distributed Snapshot
879
880The Snapshot operation behaves slightly differently when executed distributedly. First off, it requires the path
881supplied to the Snapshot call to be accessible from any worker of the cluster and from the client machine (in general
882it should be provided as an absolute path). Another important difference is that `n` separate files will be produced,
883where `n` is the number of dataset partitions. As with local RDataFrame, the result of a Snapshot on a distributed
884RDataFrame is another distributed RDataFrame on which we can define a new computation graph and run more distributed
885computations.
886
887### Distributed RunGraphs
888
889Submitting multiple distributed RDataFrame executions is supported through the RunGraphs function. Similarly to its
890local counterpart, the function expects an iterable of objects representing an RDataFrame action. Each action will be
891triggered concurrently to send multiple computation graphs to a distributed cluster at the same time:
892
893~~~{.py}
894import ROOT
895
896# Create 3 different dataframes and book an histogram on each one
897histoproxies = [
898 ROOT.RDataFrame(100, executor=SupportedExecutor(...))
899 .Define("x", "rdfentry_")
900 .Histo1D(("name", "title", 10, 0, 100), "x")
901 for _ in range(4)
902]
903
904# Execute the 3 computation graphs
905ROOT.RDF.RunGraphs(histoproxies)
906# Retrieve all the histograms in one go
907histos = [histoproxy.GetValue() for histoproxy in histoproxies]
908~~~
909
910Every distributed backend supports this feature and graphs belonging to different backends can be still triggered with
911a single call to RunGraphs (e.g. it is possible to send a Spark job and a Dask job at the same time).
912
913### Histogram models in distributed mode
914
915When calling a Histo*D operation in distributed mode, remember to pass to the function the model of the histogram to be
916computed, e.g. the axis range and the number of bins:
917
918~~~{.py}
919import ROOT
920
921if __name__ == "__main__":
922 df = ROOT.RDataFrame("mytree","myfile.root",executor=SupportedExecutor(...)).Define("x","someoperation")
923 # The model can be passed either as a tuple with the arguments in the correct order
924 df.Histo1D(("name", "title", 10, 0, 10), "x")
925 # Or by creating the specific struct
926 model = ROOT.RDF.TH1DModel("name", "title", 10, 0, 10)
927 df.Histo1D(model, "x")
928~~~
929
930Without this, two partial histograms resulting from two distributed tasks would have incompatible binning, thus leading
931to errors when merging them. Failing to pass a histogram model will raise an error on the client side, before starting
932the distributed execution.
933
934### Live visualization in distributed mode with dask
935
936The live visualization feature allows real-time data representation of plots generated during the execution
937of a distributed RDataFrame application.
938It enables visualizing intermediate results as they are computed across multiple nodes of a Dask cluster
939by creating a canvas and continuously updating it as partial results become available.
940
941The LiveVisualize() function can be imported from the Python package **ROOT.RDF.Distributed**:
942
943~~~{.py}
944import ROOT
945
946LiveVisualize = ROOT.RDF.Distributed.LiveVisualize
947~~~
948
949The function takes drawable objects (e.g. histograms) and optional callback functions as argument, it accepts 4 different input formats:
950
951- Passing a list or tuple of drawables:
952You can pass a list or tuple containing the plots you want to visualize. For example:
953
954~~~{.py}
955LiveVisualize([h_gaus, h_exp, h_random])
956~~~
957
958- Passing a list or tuple of drawables with a global callback function:
959You can also include a global callback function that will be applied to all plots. For example:
960
961~~~{.py}
962def set_fill_color(hist):
963 hist.SetFillColor("kBlue")
964
965LiveVisualize([h_gaus, h_exp, h_random], set_fill_color)
966~~~
967
968- Passing a Dictionary of drawables and callback functions:
969For more control, you can create a dictionary where keys are plots and values are corresponding (optional) callback functions. For example:
970
971~~~{.py}
972plot_callback_dict = {
973 graph: set_marker,
974 h_exp: fit_exp,
975 tprofile_2d: None
976}
977
978LiveVisualize(plot_callback_dict)
979~~~
980
981- Passing a Dictionary of drawables and callback functions with a global callback function:
982You can also combine a dictionary of plots and callbacks with a global callback function:
983
984~~~{.py}
985LiveVisualize(plot_callback_dict, write_to_tfile)
986~~~
987
988\note The allowed operations to pass to LiveVisualize are:
989 - Histo1D(), Histo2D(), Histo3D()
990 - Graph()
991 - Profile1D(), Profile2D()
992
993\warning The Live Visualization feature is only supported for the Dask backend.
994
995### Injecting C++ code and using external files into distributed RDF script
996
997Distributed RDF provides an interface for the users who want to inject the C++ code (via header files, shared libraries or declare the code directly)
998into their distributed RDF application, or their application needs to use information from external files which should be distributed
999to the workers (for example, a JSON or a txt file with necessary parameters information).
1000
1001The examples below show the usage of these interface functions: firstly, how this is done in a local Python
1002RDF application and secondly, how it is done distributedly.
1003
1004#### Include and distribute header files.
1005
1006~~~{.py}
1007# Local RDataFrame script
1008ROOT.gInterpreter.AddIncludePath("myheader.hxx")
1009df.Define(...)
1010
1011# Distributed RDF script
1012ROOT.RDF.Distributed.DistributeHeaders("myheader.hxx")
1013df.Define(...)
1014~~~
1015
1016#### Load and distribute shared libraries
1017
1018~~~{.py}
1019# Local RDataFrame script
1020ROOT.gSystem.Load("my_library.so")
1021df.Define(...)
1022
1023# Distributed RDF script
1024ROOT.RDF.Distributed.DistributeSharedLibs("my_library.so")
1025df.Define(...)
1026~~~
1027
1028#### Declare and distribute the cpp code
1029
1030The cpp code is always available to all dataframes.
1031
1032~~~{.py}
1033# Local RDataFrame script
1034ROOT.gInterpreter.Declare("my_code")
1035df.Define(...)
1036
1037# Distributed RDF script
1038ROOT.RDF.Distributed.DistributeCppCode("my_code")
1039df.Define(...)
1040~~~
1041
1042#### Distribute additional files (other than headers or shared libraries).
1043
1044~~~{.py}
1045# Local RDataFrame script is not applicable here as local RDF application can simply access the external files it needs.
1046
1047# Distributed RDF script
1048ROOT.RDF.Distributed.DistributeFiles("my_file")
1049df.Define(...)
1050~~~
1051
1052\anchor parallel-execution
1053## Performance tips and parallel execution
1054As pointed out before in this document, RDataFrame can transparently perform multi-threaded event loops to speed up
1055the execution of its actions. Users have to call ROOT::EnableImplicitMT() *before* constructing the RDataFrame
1056object to indicate that it should take advantage of a pool of worker threads. **Each worker thread processes a distinct
1057subset of entries**, and their partial results are merged before returning the final values to the user.
1058
1059By default, RDataFrame will use as many threads as the hardware supports, using up **all** the resources on
1060a machine. This might be undesirable on shared computing resources such as a batch cluster. Therefore, when running on shared computing resources, use
1061~~~{.cpp}
1062ROOT::EnableImplicitMT(numThreads)
1063~~~
1064or export an environment variable:
1065~~~{.sh}
1066export ROOT_MAX_THREADS=numThreads
1067root.exe rdfAnalysis.cxx
1068# or
1069ROOT_MAX_THREADS=4 python rdfAnalysis.py
1070~~~
1071replacing `numThreads` with the number of CPUs/slots that were allocated for this job.
1072
1073\warning There are no guarantees on the order in which threads will process the batches of
1074entries. In particular, note that this means that, for multi-thread event loops, there is no
1075guarantee on the order in which Snapshot() will _write_ entries: they could be scrambled with respect to the input
1076dataset. 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
1077runs. Likewise, Take(), AsNumpy(), ... do not preserve the original ordering.
1078
1079### Thread-safety of user-defined expressions
1080RDataFrame operations such as Histo1D() or Snapshot() are guaranteed to work correctly in multi-thread event loops.
1081User-defined expressions, such as strings or lambdas passed to Filter(), Define(), Foreach(), Reduce() or Aggregate()
1082will have to be thread-safe, i.e. it should be possible to call them concurrently from different threads.
1083
1084Note that simple Filter() and Define() transformations will inherently satisfy this requirement: Filter() / Define()
1085expressions will often be *pure* in the functional programming sense (no side-effects, no dependency on external state),
1086which eliminates all risks of race conditions.
1087
1088In order to facilitate writing of thread-safe operations, some RDataFrame features such as Foreach(), Define() or \link ROOT::RDF::RResultPtr::OnPartialResult OnPartialResult()\endlink
1089offer thread-aware counterparts (ForeachSlot(), DefineSlot(), \link ROOT::RDF::RResultPtr::OnPartialResultSlot OnPartialResultSlot()\endlink): their only difference is that they
1090will pass an extra `slot` argument (an unsigned integer) to the user-defined expression. When calling user-defined code
1091concurrently, RDataFrame guarantees that different threads will see different values of the `slot` parameter,
1092where `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.
1093In other words, within a slot, computations run sequentially, and events are processed sequentially.
1094Note that the same slot might be associated to different threads over the course of a single event loop, but two threads
1095will never receive the same slot at the same time.
1096This extra parameter might facilitate writing safe parallel code by having each thread write/modify a different
1097processing slot, e.g. a different element of a list. See [here](#generic-actions) for an example usage of ForeachSlot().
1098
1099### Parallel execution of multiple RDataFrame event loops
1100A complex analysis may require multiple separate RDataFrame computation graphs to produce all desired results. This poses the challenge that the
1101event loops of each computation graph can be parallelized, but the different loops run sequentially, one after the other.
1102On many-core architectures it might be desirable to run different event loops concurrently to improve resource usage.
1103ROOT::RDF::RunGraphs() allows running multiple RDataFrame event loops concurrently:
1104~~~{.cpp}
1105ROOT::EnableImplicitMT();
1106ROOT::RDataFrame df1("tree1", "f1.root");
1107ROOT::RDataFrame df2("tree2", "f2.root");
1108auto histo1 = df1.Histo1D("x");
1109auto histo2 = df2.Histo1D("y");
1110
1111// just accessing result pointers, the event loops of separate RDataFrames run one after the other
1112histo1->Draw(); // runs first multi-thread event loop
1113histo2->Draw(); // runs second multi-thread event loop
1114
1115// alternatively, with ROOT::RDF::RunGraphs, event loops for separate computation graphs can run concurrently
1116ROOT::RDF::RunGraphs({histo1, histo2});
1117histo1->Draw(); // results can then be used as usual
1118~~~
1119
1120### Performance considerations
1121
1122To obtain the maximum performance out of RDataFrame, make sure to avoid just-in-time compiled versions of transformations and actions if at all possible.
1123For 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.
1124Similarly, `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
1125should be preferred for performance-critical applications.
1126
1127Python applications cannot easily specify template parameters or pass C++ callables to RDataFrame.
1128See [Python interface](classROOT_1_1RDataFrame.html#python) for possible ways to speed up hot paths in this case.
1129
1130Just-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*
1131before 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.
1132
1133Also 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).
1134
1135### Memory usage
1136
1137There are two reasons why RDataFrame may consume more memory than expected.
1138
1139#### 1. Histograms in multi-threaded mode
1140In 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.
1141The 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.
1142For three-dimensional histograms, the number of clones can be reduced using ROOT::RDF::Experimental::ThreadsPerTH3().
1143~~~{.cpp}
1144#include "ROOT/RDFHelpers.hxx"
1145
1146// Make four threads share a TH3 instance:
1147ROOT::RDF::Experimental::ThreadsPerTH3(4);
1148ROOT::RDataFrame rdf(...);
1149~~~
1150
1151When 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.
1152The 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
1153slower execution and memory savings.
1154
1155#### 2. Just-in-time compilation
1156Secondly, 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:
1157
1158~~~{.cpp}
1159// assuming df1 and df2 are separate computation graphs, do:
1160auto h1 = df1.Histo1D("x");
1161auto h2 = df2.Histo1D("y");
1162h1->Draw(); // we just-in-time compile everything needed by df1 and df2 here
1163h2->Draw("SAME");
1164
1165// do not:
1166auto h1 = df1.Histo1D("x");
1167h1->Draw(); // we just-in-time compile here
1168auto h2 = df2.Histo1D("y");
1169h2->Draw("SAME"); // we just-in-time compile again here, as the second Histo1D call is new
1170~~~
1171
1172\anchor more-features
1173## More features
1174Here is a list of the most important features that have been omitted in the "Crash course" for brevity.
1175You don't need to read all these to start using RDataFrame, but they are useful to save typing time and runtime.
1176
1177\anchor systematics
1178### Systematic variations
1179
1180Starting from ROOT v6.26, RDataFrame provides a flexible syntax to define systematic variations.
1181This is done in two steps: a) register variations for one or more existing columns using Vary() and b) extract variations
1182of normal RDataFrame results using \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()". In between these steps, no other change
1183to the analysis code is required: the presence of systematic variations for certain columns is automatically propagated
1184through filters, defines and actions, and RDataFrame will take these dependencies into account when producing varied
1185results. \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" is included in header `ROOT/RDFHelpers.hxx`. The compiled C++ programs must include this header
1186explicitly, this is not required for ROOT macros.
1187
1188An example usage of Vary() and \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" in C++:
1189
1190~~~{.cpp}
1191auto nominal_hx =
1192 df.Vary("pt", "ROOT::RVecD{pt*0.9f, pt*1.1f}", {"down", "up"})
1193 .Filter("pt > pt_cut")
1194 .Define("x", someFunc, {"pt"})
1195 .Histo1D<float>("x");
1196
1197// request the generation of varied results from the nominal_hx
1198ROOT::RDF::Experimental::RResultMap<TH1D> hx = ROOT::RDF::Experimental::VariationsFor(nominal_hx);
1199
1200// the event loop runs here, upon first access to any of the results or varied results:
1201hx["nominal"].Draw(); // same effect as nominal_hx->Draw()
1202hx["pt:down"].Draw("SAME");
1203hx["pt:up"].Draw("SAME");
1204~~~
1205
1206A list of variation "tags" is passed as the last argument to Vary(). The tags give names to the varied values that are returned
1207as elements of an RVec of the appropriate C++ type. The number of variation tags must correspond to the number of elements of
1208this RVec (2 in the example above: the first element will correspond to the tag "down", the second
1209to the tag "up"). The _full_ variation name will be composed of the varied column name and the variation tags (e.g.
1210"pt:down", "pt:up" in this example). Python usage looks similar.
1211
1212Note how we use the "pt" column as usual in the Filter() and Define() calls and we simply use "x" as the value to fill
1213the resulting histogram. To produce the varied results, RDataFrame will automatically execute the Filter and Define
1214calls for each variation and fill the histogram with values and cuts that depend on the variation.
1215
1216There is no limitation to the complexity of a Vary() expression. Just like for the Define() and Filter() calls, users are
1217not limited to string expressions but they can also pass any valid C++ callable, including lambda functions and
1218complex functors. The callable can be applied to zero or more existing columns and it will always receive their
1219_nominal_ value in input.
1220
1221\note - Currently, VariationsFor() and RResultMap are in the `ROOT::RDF::Experimental` namespace, to indicate that these
1222 interfaces can still evolve and improve based on user feedback. Please send feedback on https://github.com/root-project/root.
1223
1224\note - The results of Display() cannot be varied (i.e. it is not possible to call \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()").
1225
1226See the Vary() method for more information and [this tutorial](https://root.cern/doc/master/df106__HiggsToFourLeptons_8C.html)
1227for an example usage of Vary and \ref ROOT::RDF::Experimental::VariationsFor "VariationsFor()" in the analysis.
1228
1229
1230\anchor snapshot-with-variations
1231#### Snapshot with Variations
1232To combine Vary() and Snapshot(), the \ref ROOT::RDF::RSnapshotOptions::fIncludeVariations "fIncludeVariations" flag in \ref ROOT::RDF::RSnapshotOptions
1233"RSnapshotOptions" has to be set, as demonstrated in the following example:
1234
1235~~~{.cpp}
1236 ROOT::RDF::RSnapshotOptions options;
1237 options.fIncludeVariations = true;
1238
1239 auto s = ROOT::RDataFrame(N)
1240 .Define("x", [](ULong64_t e) -> float { return 10.f * e; }, {"rdfentry_"})
1241 .Vary(
1242 "x", [](float x) { return ROOT::RVecF{x - 0.5f, x + 0.5f}; }, {"x"}, 2, "xVar")
1243 .Define("y", [](float x) -> double { return -1. * x; }, {"x"})
1244 .Filter(cuts, {"x", "y"})
1245 .Snapshot("t", filename, {"x", "y"}, options);
1246~~~
1247
1248This snapshot action will create an output file with a dataset that includes the columns requested in the Snapshot() call
1249and their corresponding variations:
1250
1251| Row | `x` | `y` | `x__xVar_0` | `y__xVar_0` | `x__xVar_1` | `y__xVar_1` | `R_rdf_mask_t_0`|
1252|------|----------------:|----------------:|----------------:|----------------:|----------------:|----------------:|----------------:|
1253| 0 | 0 | -0 | -0.5 | 0.5 | 0.5 | -0.5 | 111 |
1254| 1 | 10 | -10 | 9.5 | -9.5 | 10.5 | -10.5 | 111 |
1255| 2 | 20 | -20 | 19.5 | -19.5 | 20.5 | -20.5 | 111 |
1256| 3 | 30 | -30 | 29.5 | -29.5 | 30.5 | -30.5 | 111 |
1257| 4 | 40 | -40 | 39.5 | -39.5 | 40.5 | -40.5 | 111 |
1258| 5 | 0 | 0 | 49.5 | -49.5 | 0 | 0 | 010 |
1259| 6 |`* Not written *`| | | | | | |
1260| 7 | 0 | 0 | 0 | 0 | 70.5 | -70.5 | 100 |
1261
1262"x" and "y" are the nominal values, and their variations are snapshot as `<ColumnName>__<VariationName>_<VariationTag>`.
1263
1264When Filters are employed, some variations might not pass the selection cuts (like the rows 5 and 7 in the table above).
1265In 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
1266default-constructed object in case of classes. If none of the filters pass like in row 6, the entire event is omitted from the snapshot.
1267
1268To tell apart a genuine `0` (like `x` in row 0) from a case where nominal or variation didn't pass a selection,
1269RDataFrame writes a bitmask for each event, see last column of the table above. Every bit indicates whether its
1270associated columns are valid. The bitmask is implemented as a 64-bit `std::bitset` in memory, written to the output
1271dataset as a `std::uin64_t`. For every 64 columns, a new bitmask column is added to the output dataset.
1272
1273For each column that gets varied, the nominal and all variation columns are each assigned a bit to denote whether their
1274entries are valid. A mapping of column names to the corresponding bitmask is placed in the same file as the output
1275dataset, with a name that follows the pattern `"R_rdf_column_to_bitmask_mapping_<NAME_OF_THE_DATASET>"`. It is of type
1276`std::unordered_map<std::string, std::pair<std::string, unsigned int>>`, and maps a column name to the name of the
1277bitmask column and the index of the relevant bit. For example, in the same file as the dataset "Events" there would be
1278an object named `R_rdf_column_to_bitmask_mapping_Events`. This object for example would describe a connection such as:
1279
1280~~~
1281muon_pt --> (R_rdf_mask_Events_0, 42)
1282~~~
1283
1284which means that the validity of the entries in `muon_pt` is established by the bit `42` in the bitmask found in the
1285column `R_rdf_mask_Events_0`.
1286
1287When RDataFrame opens a file, it checks for the existence of this mapping between columns and bitmasks, and loads it automatically if found. As such,
1288RDataFrame makes the treatment of the various bitmap maskings completely transparent to the user.
1289
1290In case certain values are labeled invalid by the corresponding bit, this will result in reading a missing value. The semantics of such a scenario follow the
1291rules described in the \ref missing-values "section on dealing with missing values" and can be dealt with accordingly.
1292
1293\note Snapshot with variations is currently restricted to single-threaded TTree snapshots.
1294
1295
1296#### Varying multiple columns in lockstep
1297
1298In the following Python snippet we use the Vary() signature that allows varying multiple columns simultaneously or
1299"in lockstep":
1300
1301~~~{.python}
1302df.Vary(["pt", "eta"],
1303 "RVec<RVecF>{{pt*0.9, pt*1.1}, {eta*0.9, eta*1.1}}",
1304 variationTags=["down", "up"],
1305 variationName="ptAndEta")
1306~~~
1307
1308The expression returns an RVec of two RVecs: each inner vector contains the varied values for one column. The
1309inner vectors follow the same ordering as the column names that are passed as the first argument. Besides the variation tags, in
1310this case we also have to explicitly pass the variation name (here: "ptAndEta") as the default column name does not exist.
1311
1312The above call will produce variations "ptAndEta:down" and "ptAndEta:up".
1313
1314#### Combining multiple variations
1315
1316Even if a result depends on multiple variations, only one variation is applied at a time, i.e. there will be no result produced
1317by applying multiple systematic variations at the same time.
1318For example, in the following example snippet, the RResultMap instance `all_h` will contain keys "nominal", "pt:down",
1319"pt:up", "eta:0", "eta:1", but no "pt:up&&eta:0" or similar:
1320
1321~~~{.cpp}
1322auto df = _df.Vary("pt",
1323 "ROOT::RVecD{pt*0.9, pt*1.1}",
1324 {"down", "up"})
1325 .Vary("eta",
1326 [](float eta) { return RVecF{eta*0.9f, eta*1.1f}; },
1327 {"eta"},
1328 2);
1329
1330auto nom_h = df.Histo2D(histoModel, "pt", "eta");
1331auto all_hs = VariationsFor(nom_h);
1332all_hs.GetKeys(); // returns {"nominal", "pt:down", "pt:up", "eta:0", "eta:1"}
1333~~~
1334
1335Note how we passed the integer `2` instead of a list of variation tags to the second Vary() invocation: this is a
1336shorthand that automatically generates tags 0 to N-1 (in this case 0 and 1).
1337
1338
1339\anchor rnode
1340### RDataFrame objects as function arguments and return values
1341RDataFrame variables/nodes are relatively cheap to copy and it's possible to both pass them to (or move them into)
1342functions and to return them from functions. However, in general each dataframe node will have a different C++ type,
1343which includes all available compile-time information about what that node does. One way to cope with this complication
1344is to use template functions and/or C++14 auto return types:
1345~~~{.cpp}
1346template <typename RDF>
1347auto ApplySomeFilters(RDF df)
1348{
1349 return df.Filter("x > 0").Filter([](int y) { return y < 0; }, {"y"});
1350}
1351~~~
1352
1353A possibly simpler, C++11-compatible alternative is to take advantage of the fact that any dataframe node can be
1354converted (implicitly or via an explicit cast) to the common type ROOT::RDF::RNode:
1355~~~{.cpp}
1356// a function that conditionally adds a Range to an RDataFrame node.
1357RNode MaybeAddRange(RNode df, bool mustAddRange)
1358{
1359 return mustAddRange ? df.Range(1) : df;
1360}
1361// use as :
1362ROOT::RDataFrame df(10);
1363auto maybeRangedDF = MaybeAddRange(df, true);
1364~~~
1365
1366The conversion to ROOT::RDF::RNode is cheap, but it will introduce an extra virtual call during the RDataFrame event
1367loop (in most cases, the resulting performance impact should be negligible). Python users can perform the conversion with the helper function `ROOT.RDF.AsRNode`.
1368
1369\anchor RDFCollections
1370### Storing RDataFrame objects in collections
1371
1372ROOT::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>`:
1373
1374~~~{.cpp}
1375std::vector<ROOT::RDF::RNode> dfs;
1376dfs.emplace_back(ROOT::RDataFrame(10));
1377dfs.emplace_back(dfs[0].Define("x", "42.f"));
1378~~~
1379
1380\anchor callbacks
1381### Executing callbacks every N events
1382It's possible to schedule execution of arbitrary functions (callbacks) during the event loop.
1383Callbacks can be used e.g. to inspect partial results of the analysis while the event loop is running,
1384drawing a partially-filled histogram every time a certain number of new entries is processed, or
1385displaying a progress bar while the event loop runs.
1386
1387For example one can draw an up-to-date version of a result histogram every 100 entries like this:
1388~~~{.cpp}
1389auto h = df.Histo1D("x");
1390TCanvas c("c","x hist");
1391h.OnPartialResult(100, [&c](TH1D &h_) { c.cd(); h_.Draw(); c.Update(); });
1392// event loop runs here, this final `Draw` is executed after the event loop is finished
1393h->Draw();
1394~~~
1395
1396Callbacks are registered to a ROOT::RDF::RResultPtr and must be callables that takes a reference to the result type as argument
1397and return nothing. RDataFrame will invoke registered callbacks passing partial action results as arguments to them
1398(e.g. a histogram filled with a part of the selected events).
1399
1400Read more on ROOT::RDF::RResultPtr::OnPartialResult() and ROOT::RDF::RResultPtr::OnPartialResultSlot().
1401
1402\anchor default-branches
1403### Default column lists
1404When constructing an RDataFrame object, it is possible to specify a **default column list** for your analysis, in the
1405usual form of a list of strings representing branch/column names. The default column list will be used as a fallback
1406whenever a list specific to the transformation/action is not present. RDataFrame will take as many of these columns as
1407needed, ignoring trailing extra names if present.
1408~~~{.cpp}
1409// use "b1" and "b2" as default columns
1410ROOT::RDataFrame d1("myTree", "file.root", {"b1","b2"});
1411auto h = d1.Filter([](int b1, int b2) { return b1 > b2; }) // will act on "b1" and "b2"
1412 .Histo1D(); // will act on "b1"
1413
1414// just one default column this time
1415ROOT::RDataFrame d2("myTree", "file.root", {"b1"});
1416auto d2f = d2.Filter([](double b2) { return b2 > 0; }, {"b2"}) // we can still specify non-default column lists
1417auto min = d2f.Min(); // returns the minimum value of "b1" for the filtered entries
1418auto vals = d2f.Take<double>(); // return the values for all entries passing the selection as a vector
1419~~~
1420
1421\anchor helper-cols
1422### Special helper columns: rdfentry_ and rdfslot_
1423Every instance of RDataFrame is created with two special columns called `rdfentry_` and `rdfslot_`. The `rdfentry_`
1424column is of type `ULong64_t` and it holds the current entry number while `rdfslot_` is an `unsigned int`
1425holding the index of the current data processing slot.
1426For backwards compatibility reasons, the names `tdfentry_` and `tdfslot_` are also accepted.
1427These columns are ignored by operations such as [Cache](classROOT_1_1RDF_1_1RInterface.html#aaaa0a7bb8eb21315d8daa08c3e25f6c9)
1428or [Snapshot](classROOT_1_1RDF_1_1RInterface.html#a233b7723e498967f4340705d2c4db7f8).
1429
1430\warning Note that in multi-thread event loops the values of `rdfentry_` _do not_ correspond to what would be the entry numbers
1431of a TChain constructed over the same set of ROOT files, as the entries are processed in an unspecified order.
1432
1433\anchor jitting
1434### Just-in-time compilation: column type inference and explicit declaration of column types
1435C++ is a statically typed language: all types must be known at compile-time. This includes the types of the TTree
1436branches we want to work on. For filters, defined columns and some of the actions, **column types are deduced from the
1437signature** of the relevant filter function/temporary column expression/action function:
1438~~~{.cpp}
1439// here b1 is deduced to be `int` and b2 to be `double`
1440df.Filter([](int x, double y) { return x > 0 && y < 0.; }, {"b1", "b2"});
1441~~~
1442If we specify an incorrect type for one of the columns, an exception with an informative message will be thrown at
1443runtime, when the column value is actually read from the dataset: RDataFrame detects type mismatches. The same would
1444happen if we swapped the order of "b1" and "b2" in the column list passed to Filter().
1445
1446Certain actions, on the other hand, do not take a function as argument (e.g. Histo1D()), so we cannot deduce the type of
1447the column at compile-time. In this case **RDataFrame infers the type of the column** from the TTree itself. This
1448is why we never needed to specify the column types for all actions in the above snippets.
1449
1450When the column type is not a common one such as `int`, `double`, `char` or `float` it is nonetheless good practice to
1451specify it as a template parameter to the action itself, like this:
1452~~~{.cpp}
1453df.Histo1D("b1"); // OK, the type of "b1" is deduced at runtime
1454df.Min<MyNumber_t>("myObject"); // OK, "myObject" is deduced to be of type `MyNumber_t`
1455~~~
1456
1457Deducing types at runtime requires the just-in-time compilation of the relevant actions, which has a small runtime
1458overhead, so specifying the type of the columns as template parameters to the action is good practice when performance is a goal.
1459
1460When 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`:
1461
1462~~~{.cpp}
1463// this throws an error (note the typo)
1464df.Define("x", "0").Filter("x = 0");
1465~~~
1466
1467\anchor generic-actions
1468### User-defined custom actions
1469RDataFrame strives to offer a comprehensive set of standard actions that can be performed on each event. At the same
1470time, it allows users to inject their own action code to perform arbitrarily complex data reductions.
1471
1472#### Implementing custom actions with Book()
1473
1474Through the Book() method, users can implement a custom action and have access to the same features
1475that built-in RDataFrame actions have, e.g. hooks to events related to the start, end and execution of the
1476event loop, or the possibility to return a lazy RResultPtr to an arbitrary type of result:
1477
1478~~~{.cpp}
1479#include <ROOT/RDataFrame.hxx>
1480#include <memory>
1481
1482class MyCounter : public ROOT::Detail::RDF::RActionImpl<MyCounter> {
1483 std::shared_ptr<int> fFinalResult = std::make_shared<int>(0);
1484 std::vector<int> fPerThreadResults;
1485
1486public:
1487 // We use a public type alias to advertise the type of the result of this action
1488 using Result_t = int;
1489
1490 MyCounter(unsigned int nSlots) : fPerThreadResults(nSlots) {}
1491
1492 // Called before the event loop to retrieve the address of the result that will be filled/generated.
1493 std::shared_ptr<int> GetResultPtr() const { return fFinalResult; }
1494
1495 // Called at the beginning of the event loop.
1496 void Initialize() {}
1497
1498 // Called at the beginning of each processing task.
1499 void InitTask(TTreeReader *, int) {}
1500
1501 /// Called at every entry.
1502 void Exec(unsigned int slot)
1503 {
1504 fPerThreadResults[slot]++;
1505 }
1506
1507 // Called at the end of the event loop.
1508 void Finalize()
1509 {
1510 *fFinalResult = std::accumulate(fPerThreadResults.begin(), fPerThreadResults.end(), 0);
1511 }
1512
1513 // Called by RDataFrame to retrieve the name of this action.
1514 std::string GetActionName() const { return "MyCounter"; }
1515};
1516
1517int main() {
1518 ROOT::RDataFrame df(10);
1519 ROOT::RDF::RResultPtr<int> resultPtr = df.Book<>(MyCounter{df.GetNSlots()}, {});
1520 // The GetValue call triggers the event loop
1521 std::cout << "Number of processed entries: " << resultPtr.GetValue() << std::endl;
1522}
1523~~~
1524
1525See the Book() method for more information and [this tutorial](https://root.cern/doc/master/df018__customActions_8C.html)
1526for a more complete example.
1527
1528#### Injecting arbitrary code in the event loop with Foreach() and ForeachSlot()
1529
1530Foreach() takes a callable (lambda expression, free function, functor...) and a list of columns and
1531executes the callable on the values of those columns for each event that passes all upstream selections.
1532It can be used to perform actions that are not already available in the interface. For example, the following snippet
1533evaluates the root mean square of column "x":
1534~~~{.cpp}
1535// Single-thread evaluation of RMS of column "x" using Foreach
1536double sumSq = 0.;
1537unsigned int n = 0;
1538df.Foreach([&sumSq, &n](double x) { ++n; sumSq += x*x; }, {"x"});
1539std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1540~~~
1541In multi-thread runs, users are responsible for the thread-safety of the expression passed to Foreach():
1542thread will execute the expression concurrently.
1543The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of `rms`; but
1544this is probably too much head-scratch for such a simple operation.
1545
1546ForeachSlot() can help in this situation. It is an alternative version of Foreach() for which the function takes an
1547additional "processing slot" parameter besides the columns it should be applied to. RDataFrame
1548guarantees that ForeachSlot() will invoke the user expression with different `slot` parameters for different concurrent
1549executions (see [Special helper columns: rdfentry_ and rdfslot_](\ref helper-cols) for more information on the slot parameter).
1550We can take advantage of ForeachSlot() to evaluate a thread-safe root mean square of column "x":
1551~~~{.cpp}
1552// Thread-safe evaluation of RMS of column "x" using ForeachSlot
1553ROOT::EnableImplicitMT();
1554const unsigned int nSlots = df.GetNSlots();
1555std::vector<double> sumSqs(nSlots, 0.);
1556std::vector<unsigned int> ns(nSlots, 0);
1557
1558df.ForeachSlot([&sumSqs, &ns](unsigned int slot, double x) { sumSqs[slot] += x*x; ns[slot] += 1; }, {"x"});
1559double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
1560unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
1561std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1562~~~
1563Notice how we created one `double` variable for each processing slot and later merged their results via `std::accumulate`.
1564
1565
1566\anchor friends
1567### Dataset joins with friend trees
1568
1569Vertically concatenating multiple trees that have the same columns (creating a logical dataset with the same columns and
1570more rows) is trivial in RDataFrame: just pass the tree name and a list of file names to RDataFrame's constructor, or create a TChain
1571out of the desired trees and pass that to RDataFrame.
1572
1573Horizontal concatenations of trees or chains (creating a logical dataset with the same number of rows and the union of the
1574columns of multiple trees) leverages TTree's "friend" mechanism.
1575
1576Simple joins of trees that do not have the same number of rows are also possible with indexed friend trees (see below).
1577
1578To use friend trees in RDataFrame, set up trees with the appropriate relationships and then instantiate an RDataFrame
1579with the main tree:
1580
1581~~~{.cpp}
1582TTree main([...]);
1583TTree friend([...]);
1584main.AddFriend(&friend, "myFriend");
1585
1586RDataFrame df(main);
1587auto df2 = df.Filter("myFriend.MyCol == 42");
1588~~~
1589
1590The same applies for TChains. Columns coming from the friend trees can be referred to by their full name, like in the example above,
1591or the friend tree name can be omitted in case the column name is not ambiguous (e.g. "MyCol" could be used instead of
1592"myFriend.MyCol" in the example above if there is no column "MyCol" in the main tree).
1593
1594\note A common source of confusion is that trees that are written out from a multi-thread Snapshot() call will have their
1595 entries (block-wise) shuffled with respect to the original tree. Such trees cannot be used as friends of the original
1596 one: rows will be mismatched.
1597
1598Indexed friend trees provide a way to perform simple joins of multiple trees over a common column.
1599When a certain entry in the main tree (or chain) is loaded, the friend trees (or chains) will then load an entry where the
1600"index" columns have a value identical to the one in the main one. For example, in Python:
1601
1602~~~{.py}
1603main_tree = ...
1604aux_tree = ...
1605
1606# If a friend tree has an index on `commonColumn`, when the main tree loads
1607# a given row, it also loads the row of the friend tree that has the same
1608# value of `commonColumn`
1609aux_tree.BuildIndex("commonColumn")
1610
1611mainTree.AddFriend(aux_tree)
1612
1613df = ROOT.RDataFrame(mainTree)
1614~~~
1615
1616RDataFrame supports indexed friend TTrees from ROOT v6.24 in single-thread mode and from v6.28/02 in multi-thread mode.
1617
1618\anchor other-file-formats
1619### Reading data formats other than ROOT trees
1620RDataFrame can be interfaced with RDataSources. The ROOT::RDF::RDataSource interface defines an API that RDataFrame can use to read arbitrary columnar data formats.
1621
1622RDataFrame calls into concrete RDataSource implementations to retrieve information about the data, retrieve (thread-local) readers or "cursors" for selected columns
1623and to advance the readers to the desired data entry.
1624Some predefined RDataSources are natively provided by ROOT such as the ROOT::RDF::RCsvDS which allows to read comma separated files:
1625~~~{.cpp}
1626auto tdf = ROOT::RDF::FromCSV("MuRun2010B.csv");
1627auto filteredEvents =
1628 tdf.Filter("Q1 * Q2 == -1")
1629 .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
1630auto h = filteredEvents.Histo1D("m");
1631h->Draw();
1632~~~
1633
1634See also FromNumpy (Python-only), FromRNTuple(), FromArrow(), FromSqlite().
1635
1636\anchor callgraphs
1637### Computation graphs (storing and reusing sets of transformations)
1638
1639As 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
1640several paths of filtering/creation of columns are executed simultaneously, and finally aggregated results are produced.
1641
1642RDataFrame detects when several actions use the same filter or the same defined column, and **only evaluates each
1643filter or defined column once per event**, regardless of how many times that result is used down the computation graph.
1644Objects read from each column are **built once and never copied**, for maximum efficiency.
1645When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated,
1646so it might be advisable to put the strictest filters first in the graph.
1647
1648\anchor representgraph
1649### Visualizing the computation graph
1650It 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
1651or in a file.
1652
1653Invoking the function ROOT::RDF::SaveGraph() on any node that is not the head node, the computation graph of the branch
1654the node belongs to is printed. By using the head node, the entire computation graph is printed.
1655
1656Following there is an example of usage:
1657~~~{.cpp}
1658// First, a sample computational graph is built
1659ROOT::RDataFrame df("tree", "f.root");
1660
1661auto df2 = df.Define("x", []() { return 1; })
1662 .Filter("col0 % 1 == col0")
1663 .Filter([](int b1) { return b1 <2; }, {"cut1"})
1664 .Define("y", []() { return 1; });
1665
1666auto count = df2.Count();
1667
1668// Prints the graph to the rd1.dot file in the current directory
1669ROOT::RDF::SaveGraph(df, "./mydot.dot");
1670// Prints the graph to standard output
1671ROOT::RDF::SaveGraph(df);
1672~~~
1673
1674The 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:
1675~~~{.sh}
1676$ dot -Tpng computation_graph.dot -ocomputation_graph.png
1677~~~
1678
1679\image html RDF_Graph2.png
1680
1681\anchor rdf-logging
1682### Activating RDataFrame execution logs
1683
1684RDataFrame has experimental support for verbose logging of the event loop runtimes and other interesting related information. It is activated as follows:
1685~~~{.cpp}
1686#include <ROOT/RLogger.hxx>
1687
1688// this increases RDF's verbosity level as long as the `verbosity` variable is in scope
1689auto verbosity = ROOT::RLogScopedVerbosity(ROOT::Detail::RDF::RDFLogChannel(), ROOT::ELogLevel::kInfo);
1690~~~
1691
1692or in Python:
1693~~~{.python}
1694import ROOT
1695
1696verbosity = ROOT.RLogScopedVerbosity(ROOT.Detail.RDF.RDFLogChannel(), ROOT.ELogLevel.kInfo)
1697~~~
1698
1699More information (e.g. start and end of each multi-thread task) is printed using `ELogLevel.kDebug` and even more
1700(e.g. a full dump of the generated code that RDataFrame just-in-time-compiles) using `ELogLevel.kDebug+10`.
1701
1702\anchor rdf-from-spec
1703### Creating an RDataFrame from a dataset specification file
1704
1705RDataFrame can be created using a dataset specification JSON file:
1706
1707~~~{.python}
1708import ROOT
1709
1710df = ROOT.RDF.Experimental.FromSpec("spec.json")
1711~~~
1712
1713The input dataset specification JSON file needs to be provided by the user and it describes all necessary samples and
1714their associated metadata information. The main required key is the "samples" (at least one sample is needed) and the
1715required sub-keys for each sample are "trees" and "files". Additionally, one can specify a metadata dictionary for each
1716sample in the "metadata" key.
1717
1718A simple example for the formatting of the specification in the JSON file is the following:
1719
1720~~~{.cpp}
1721{
1722 "samples": {
1723 "sampleA": {
1724 "trees": ["tree1", "tree2"],
1725 "files": ["file1.root", "file2.root"],
1726 "metadata": {
1727 "lumi": 10000.0,
1728 "xsec": 1.0,
1729 "sample_category" = "data"
1730 }
1731 },
1732 "sampleB": {
1733 "trees": ["tree3", "tree4"],
1734 "files": ["file3.root", "file4.root"],
1735 "metadata": {
1736 "lumi": 0.5,
1737 "xsec": 1.5,
1738 "sample_category" = "MC_background"
1739 }
1740 }
1741 }
1742}
1743~~~
1744
1745The metadata information from the specification file can be then accessed using the DefinePerSample function.
1746For example, to access luminosity information (stored as a double):
1747
1748~~~{.python}
1749df.DefinePerSample("lumi", 'rdfsampleinfo_.GetD("lumi")')
1750~~~
1751
1752or sample_category information (stored as a string):
1753
1754~~~{.python}
1755df.DefinePerSample("sample_category", 'rdfsampleinfo_.GetS("sample_category")')
1756~~~
1757
1758or directly the filename:
1759
1760~~~{.python}
1761df.DefinePerSample("name", "rdfsampleinfo_.GetSampleName()")
1762~~~
1763
1764An example implementation of the "FromSpec" method is available in tutorial: df106_HiggstoFourLeptons.py, which also
1765provides a corresponding exemplary JSON file for the dataset specification.
1766
1767\anchor progressbar
1768### Adding a progress bar
1769
1770A progress bar showing the processed event statistics can be added to any RDataFrame program.
1771The event statistics include elapsed time, currently processed file, currently processed events, the rate of event processing
1772and an estimated remaining time (per file being processed). It is recorded and printed in the terminal every m events and every
1773n seconds (by default m = 1000 and n = 1). The ProgressBar can be also added when the multithread (MT) mode is enabled.
1774
1775ProgressBar is added after creating the dataframe object (df):
1776~~~{.cpp}
1777ROOT::RDataFrame df("tree", "file.root");
1778ROOT::RDF::Experimental::AddProgressBar(df);
1779~~~
1780
1781Alternatively, RDataFrame can be cast to an RNode first, giving the user more flexibility
1782For example, it can be called at any computational node, such as Filter or Define, not only the head node,
1783with no change to the ProgressBar function itself (please see the [Python interface](classROOT_1_1RDataFrame.html#python)
1784section for appropriate usage in Python):
1785~~~{.cpp}
1786ROOT::RDataFrame df("tree", "file.root");
1787auto df_1 = ROOT::RDF::RNode(df.Filter("x>1"));
1788ROOT::RDF::Experimental::AddProgressBar(df_1);
1789~~~
1790Examples 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).
1791
1792\anchor missing-values
1793### Working with missing values in the dataset
1794
1795In certain situations a dataset might be missing one or more values at one or
1796more of its entries. For example:
1797
1798- If the dataset is composed of multiple files and one or more files is
1799 missing one or more columns required by the analysis.
1800- When joining different datasets horizontally according to some index value
1801 (e.g. the event number), if the index does not find a match in one or more
1802 other datasets for a certain entry.
1803- If, for a certain event, a column is invalid because it results from a Snapshot
1804 with systematic variations, and that variation didn't pass its filters. For
1805 more details, see \ref snapshot-with-variations.
1806
1807For example, suppose that column "y" does not have a value for entry 42:
1808
1809\code{.unparsed}
1810+-------+---+---+
1811| Entry | x | y |
1812+-------+---+---+
1813| 42 | 1 | |
1814+-------+---+---+
1815\endcode
1816
1817If the RDataFrame application reads that column, for example if a Take() action
1818was requested, the default behaviour is to throw an exception indicating
1819that that column is missing an entry.
1820
1821The following paragraphs discuss the functionalities provided by RDataFrame to
1822work with missing values in the dataset.
1823
1824#### FilterAvailable and FilterMissing
1825
1826FilterAvailable and FilterMissing are specialized RDataFrame Filter operations.
1827They take as input argument the name of a column of the dataset to watch for
1828missing values. Like Filter, they will either keep or discard an entire entry
1829based on whether a condition returns true or false. Specifically:
1830
1831- FilterAvailable: the condition is whether the value of the column is present.
1832 If so, the entry is kept. Otherwise if the value is missing the entry is
1833 discarded.
1834- FilterMissing: the condition is whether the value of the column is missing. If
1835 so, the entry is kept. Otherwise if the value is present the entry is
1836 discarded.
1837
1838\code{.py}
1839df = ROOT.RDataFrame(dataset)
1840
1841# Anytime an entry from "col" is missing, the entire entry will be filtered out
1842df_available = df.FilterAvailable("col")
1843df_available = df_available.Define("twice", "col * 2")
1844
1845# Conversely, if we want to select the entries for which the column has missing
1846# values, we do the following
1847df_missingcol = df.FilterMissing("col")
1848# Following operations in the same branch of the computation graph clearly
1849# cannot access that same column, since there would be no value to read
1850df_missingcol = df_missingcol.Define("observable", "othercolumn * 2")
1851\endcode
1852
1853\code{.cpp}
1854ROOT::RDataFrame df{dataset};
1855
1856// Anytime an entry from "col" is missing, the entire entry will be filtered out
1857auto df_available = df.FilterAvailable("col");
1858auto df_twicecol = df_available.Define("twice", "col * 2");
1859
1860// Conversely, if we want to select the entries for which the column has missing
1861// values, we do the following
1862auto df_missingcol = df.FilterMissing("col");
1863// Following operations in the same branch of the computation graph clearly
1864// cannot access that same column, since there would be no value to read
1865auto df_observable = df_missingcol.Define("observable", "othercolumn * 2");
1866\endcode
1867
1868#### DefaultValueFor
1869
1870DefaultValueFor creates a node of the computation graph which just forwards the
1871values of the columns necessary for other downstream nodes, when they are
1872available. In case a value of the input column passed to this function is not
1873available, the node will provide the default value passed to this function call
1874instead. Example:
1875
1876\code{.py}
1877df = ROOT.RDataFrame(dataset)
1878# Anytime an entry from "col" is missing, the value will be the default one
1879default_value = ... # Some sensible default value here
1880df = df.DefaultValueFor("col", default_value)
1881df = df.Define("twice", "col * 2")
1882\endcode
1883
1884\code{.cpp}
1885ROOT::RDataFrame df{dataset};
1886// Anytime an entry from "col" is missing, the value will be the default one
1887constexpr auto default_value = ... // Some sensible default value here
1888auto df_default = df.DefaultValueFor("col", default_value);
1889auto df_col = df_default.Define("twice", "col * 2");
1890\endcode
1891
1892#### Mixing different strategies to work with missing values in the same RDataFrame
1893
1894All the operations presented above only act on the particular branch of the
1895computation graph where they are called, so that different results can be
1896obtained by mixing and matching the filtering or providing a default value
1897strategies:
1898
1899\code{.py}
1900df = ROOT.RDataFrame(dataset)
1901# Anytime an entry from "col" is missing, the value will be the default one
1902default_value = ... # Some sensible default value here
1903df_default = df.DefaultValueFor("col", default_value).Define("twice", "col * 2")
1904df_filtered = df.FilterAvailable("col").Define("twice", "col * 2")
1905
1906# Same number of total entries as the input dataset, with defaulted values
1907df_default.Display(["twice"]).Print()
1908# Only keep the entries where "col" has values
1909df_filtered.Display(["twice"]).Print()
1910\endcode
1911
1912\code{.cpp}
1913ROOT::RDataFrame df{dataset};
1914
1915// Anytime an entry from "col" is missing, the value will be the default one
1916constexpr auto default_value = ... // Some sensible default value here
1917auto df_default = df.DefaultValueFor("col", default_value).Define("twice", "col * 2");
1918auto df_filtered = df.FilterAvailable("col").Define("twice", "col * 2");
1919
1920// Same number of total entries as the input dataset, with defaulted values
1921df_default.Display({"twice"})->Print();
1922// Only keep the entries where "col" has values
1923df_filtered.Display({"twice"})->Print();
1924\endcode
1925
1926#### Further considerations
1927
1928Note that working with missing values is currently supported with a TTree-based
1929data source. Support of this functionality for other data sources may come in
1930the future.
1931
1932\anchor special-values
1933### Dealing with NaN or Inf values in the dataset
1934
1935RDataFrame does not treat NaNs or infinities beyond what the floating-point standards require, i.e. they will
1936propagate to the final result.
1937Non-finite numbers can be suppressed using Filter(), e.g.:
1938
1939\code{.py}
1940df.Filter("std::isfinite(x)").Mean("x")
1941\endcode
1942
1943*/
1944// clang-format on
1945
1946namespace ROOT {
1947
1949using ColumnNamesPtr_t = std::shared_ptr<const ColumnNames_t>;
1950
1951////////////////////////////////////////////////////////////////////////////
1952/// \brief Build the dataframe.
1953/// \param[in] treeName Name of the tree contained in the directory
1954/// \param[in] dirPtr TDirectory where the tree is stored, e.g. a TFile.
1955/// \param[in] defaultColumns Collection of default columns.
1956///
1957/// The default columns are looked at in case no column is specified in the
1958/// booking of actions or transformations.
1959/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
1960RDataFrame::RDataFrame(std::string_view treeName, TDirectory *dirPtr, const ColumnNames_t &defaultColumns)
1961 : RInterface(std::make_shared<RDFDetail::RLoopManager>(
1962 std::make_unique<ROOT::Internal::RDF::RTTreeDS>(treeName, dirPtr), defaultColumns))
1963{
1964}
1965
1966////////////////////////////////////////////////////////////////////////////
1967/// \brief Build the dataframe.
1968/// \param[in] treeName Name of the tree contained in the directory
1969/// \param[in] fileNameGlob TDirectory where the tree is stored, e.g. a TFile.
1970/// \param[in] defaultColumns Collection of default columns.
1971///
1972/// The filename glob supports the same type of expressions as TChain::Add(), and it is passed as-is to TChain's
1973/// constructor.
1974///
1975/// The default columns are looked at in case no column is specified in the
1976/// booking of actions or transformations.
1977/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
1978RDataFrame::RDataFrame(std::string_view treeName, std::string_view fileNameGlob, const ColumnNames_t &defaultColumns)
1979 : RInterface(ROOT::Detail::RDF::CreateLMFromFile(treeName, fileNameGlob, defaultColumns))
1980{
1981}
1982
1983////////////////////////////////////////////////////////////////////////////
1984/// \brief Build the dataframe.
1985/// \param[in] datasetName Name of the dataset contained in the directory
1986/// \param[in] fileNameGlobs Collection of file names of filename globs
1987/// \param[in] defaultColumns Collection of default columns.
1988///
1989/// The filename globs support the same type of expressions as TChain::Add(), and each glob is passed as-is
1990/// to TChain's constructor.
1991///
1992/// The default columns are looked at in case no column is specified in the booking of actions or transformations.
1993/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
1994RDataFrame::RDataFrame(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1996 : RInterface(ROOT::Detail::RDF::CreateLMFromFile(datasetName, fileNameGlobs, defaultColumns))
1997{
1998}
1999
2000////////////////////////////////////////////////////////////////////////////
2001/// \brief Build the dataframe.
2002/// \param[in] tree The tree or chain to be studied.
2003/// \param[in] defaultColumns Collection of default column names to fall back to when none is specified.
2004///
2005/// The default columns are looked at in case no column is specified in the
2006/// booking of actions or transformations.
2007/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
2012
2013//////////////////////////////////////////////////////////////////////////
2014/// \brief Build a dataframe that generates numEntries entries.
2015/// \param[in] numEntries The number of entries to generate.
2016///
2017/// An empty-source dataframe constructed with a number of entries will
2018/// generate those entries on the fly when some action is triggered,
2019/// and it will do so for all the previously-defined columns.
2020/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
2022 : RInterface(std::make_shared<RDFDetail::RLoopManager>(numEntries))
2023
2024{
2025}
2026
2027//////////////////////////////////////////////////////////////////////////
2028/// \brief Build dataframe associated to data source.
2029/// \param[in] ds The data source object.
2030/// \param[in] defaultColumns Collection of default column names to fall back to when none is specified.
2031///
2032/// A dataframe associated to a data source will query it to access column values.
2033/// \note see ROOT::RDF::RInterface for the documentation of the methods available.
2034RDataFrame::RDataFrame(std::unique_ptr<ROOT::RDF::RDataSource> ds, const ColumnNames_t &defaultColumns)
2035 : RInterface(std::make_shared<RDFDetail::RLoopManager>(std::move(ds), defaultColumns))
2036{
2037}
2038
2039//////////////////////////////////////////////////////////////////////////
2040/// \brief Build dataframe from an RDatasetSpec object.
2041/// \param[in] spec The dataset specification object.
2042///
2043/// A dataset specification includes trees and file names,
2044/// as well as an optional friend list and/or entry range.
2045///
2046/// ### Example usage from Python:
2047/// ~~~{.py}
2048/// spec = (
2049/// ROOT.RDF.Experimental.RDatasetSpec()
2050/// .AddSample(("data", "tree", "file.root"))
2051/// .WithGlobalFriends("friendTree", "friend.root", "alias")
2052/// .WithGlobalRange((100, 200))
2053/// )
2054/// df = ROOT.RDataFrame(spec)
2055/// ~~~
2056///
2057/// See also ROOT::RDataFrame::FromSpec().
2062
2064{
2065 // If any node of the computation graph associated with this RDataFrame
2066 // declared code to jit, we need to make sure the compilation actually
2067 // happens. For example, a jitted Define could have been booked but
2068 // if the computation graph is not actually run then the code of the
2069 // Define node is not jitted. This in turn would cause memory leaks.
2070 // See https://github.com/root-project/root/issues/15399
2071 fLoopManager->Jit();
2072}
2073
2074namespace RDF {
2075namespace Experimental {
2076
2077////////////////////////////////////////////////////////////////////////////
2078/// \brief Create the RDataFrame from the dataset specification file.
2079/// \param[in] jsonFile Path to the dataset specification JSON file.
2080///
2081/// The input dataset specification JSON file must include a number of keys that
2082/// describe all the necessary samples and their associated metadata information.
2083/// The main key, "samples", is required and at least one sample is needed. Each
2084/// sample must have at least one key "trees" and at least one key "files" from
2085/// which the data is read. Optionally, one or more metadata information can be
2086/// added, as well as the friend list information.
2087///
2088/// ### Example specification file JSON:
2089/// The following is an example of the dataset specification JSON file formatting:
2090///~~~{.cpp}
2091/// {
2092/// "samples": {
2093/// "sampleA": {
2094/// "trees": ["tree1", "tree2"],
2095/// "files": ["file1.root", "file2.root"],
2096/// "metadata": {"lumi": 1.0, }
2097/// },
2098/// "sampleB": {
2099/// "trees": ["tree3", "tree4"],
2100/// "files": ["file3.root", "file4.root"],
2101/// "metadata": {"lumi": 0.5, }
2102/// },
2103/// ...
2104/// },
2105/// }
2106///~~~
2111
2112} // namespace Experimental
2113} // namespace RDF
2114
2115} // namespace ROOT
2116
2117namespace cling {
2118//////////////////////////////////////////////////////////////////////////
2119/// Print an RDataFrame at the prompt
2120std::string printValue(ROOT::RDataFrame *df)
2121{
2122 // The loop manager is never null, except when its construction failed.
2123 // This can happen e.g. if the constructor of RLoopManager that expects
2124 // a file name is used and that file doesn't exist. This point is usually
2125 // not even reached in that situation, since the exception thrown by the
2126 // constructor will also stop execution of the program. But it can still
2127 // be reached at the prompt, if the user tries to print the RDataFrame
2128 // variable after an incomplete initialization.
2129 auto *lm = df->GetLoopManager();
2130 if (!lm) {
2131 throw std::runtime_error("Cannot print information about this RDataFrame, "
2132 "it was not properly created. It must be discarded.");
2133 }
2134 auto defCols = lm->GetDefaultColumnNames();
2135
2136 std::ostringstream ret;
2137 if (auto ds = df->GetDataSource()) {
2138 ret << "A data frame associated to the data source \"" << cling::printValue(ds) << "\"";
2139 } else {
2140 ret << "An empty data frame that will create " << lm->GetNEmptyEntries() << " entries\n";
2141 }
2142
2143 return ret.str();
2144}
2145} // 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