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