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