Logo ROOT   6.10/09
Reference Guide
TDataFrame.cxx
Go to the documentation of this file.
1 // Author: Enrico Guiraud, Danilo Piparo CERN 12/2016
2 
3 /*************************************************************************
4  * Copyright (C) 1995-2016, 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 "ROOT/TDataFrame.hxx"
12 using namespace ROOT::Experimental;
13 
14 /**
15 * \class ROOT::Experimental::TDataFrame
16 * \ingroup dataframe
17 * \brief ROOT's TDataFrame offers a high level interface for analyses of data stored in TTrees.
18 
19 In addition, multi-threading and other low-level optimisations allow users to exploit all the resources available
20 on their machines completely transparently.<br>
21 Skip to the [class reference](#reference) or keep reading for the user guide.
22 
23 In a nutshell:
24 ~~~{.cpp}
25 ROOT::EnableImplicitMT(); // Tell ROOT you want to go parallel
26 ROOT::Experimental::TDataFrame d("myTree", file); // Interface to TTree and TChain
27 auto myHisto = d.Histo1D("Branch_A"); // This happens in parallel!
28 myHisto->Draw();
29 ~~~
30 
31 Calculations are expressed in terms of a type-safe *functional chain of actions and transformations*, `TDataFrame` takes
32 care of their execution. The implementation automatically puts in place several low level optimisations such as
33 multi-thread parallelisation and caching.
34 The namespace containing the TDataFrame is ROOT::Experimental. This signals the fact that the interfaces may evolve in
35 time.
36 
37 \htmlonly
38 <a href="https://doi.org/10.5281/zenodo.260230"><img src="https://zenodo.org/badge/DOI/10.5281/zenodo.260230.svg"
39 alt="DOI"></a>
40 \endhtmlonly
41 
42 ## Table of Contents
43 - [Introduction](#introduction)
44 - [Crash course](#crash-course)
45 - [More features](#more-features)
46 - [Transformations](#transformations)
47 - [Actions](#actions)
48 - [Parallel execution](#parallel-execution)
49 
50 ## <a name="introduction"></a>Introduction
51 A pipeline of operations is described to be performed on the data, the framework takes care
52 of the management of the loop over entries as well as low-level details such as I/O and parallelisation.
53 `TDataFrame` provides an interface to perform most common operations required by ROOT analyses;
54 at the same time, the users are not limited to those
55 common operations: building blocks to trigger custom calculations are available too.
56 
57 `TDataFrame` is built with a *modular* and *flexible* workflow in mind, summarised as follows:
58 
59 1. **build a data-frame** object by specifying your data-set
60 2. **apply a series of transformations** to your data
61  1. **filter** (e.g. apply some cuts) or
62  2. create a **temporary column** (e.g. the result of an expensive computation on branches, or an alias for a
63 branch)
64 3. **apply actions** to the transformed data to produce results (e.g. fill a histogram)
65 4.
66 <table>
67 <tr>
68  <td>
69  <b>TTreeReader</b>
70  </td>
71  <td>
72  <b>ROOT::Experimental::TDataFrame</b>
73  </td>
74 </tr>
75 <tr>
76  <td>
77 ~~~{.cpp}
78 TTreeReader reader("myTree", file);
79 TTreeReaderValue<A_t> a(reader, "A");
80 TTreeReaderValue<B_t> b(reader, "B");
81 TTreeReaderValue<C_t> c(reader, "C");
82 while(reader.Next()) {
83  if(IsGoodEvent(a, b, c))
84  DoStuff(a, b, c);
85 }
86 ~~~
87  </td>
88  <td>
89 ~~~{.cpp}
90 ROOT::Experimental::TDataFrame d("myTree", file, {"A", "B", "C"});
91 d.Filter(IsGoodEvent).Foreach(DoStuff);
92 ~~~
93  </td>
94 </tr>
95 <tr>
96  <td>
97  <b>TTree::Draw</b>
98  </td>
99  <td>
100  <b>ROOT::Experimental::TDataFrame</b>
101  </td>
102 </tr>
103 <tr>
104  <td>
105 ~~~{.cpp}
106 TTree *t = static_cast<TTree*>(
107  file->Get("myTree")
108 );
109 t->Draw("var", "var > 2");
110 ~~~
111  </td>
112  <td>
113 ~~~{.cpp}
114 ROOT::Experimental::TDataFrame d("myTree", file, "var");
115 d.Filter([](int v) { return v > 2; }).Histo1D();
116 ~~~
117  </td>
118 </tr>
119 </table>
120 
121 Keep reading to follow a five-minute [crash course](#crash-course) to `TDataFrame`, or jump to an overview of useful
122 [features](#more-features), or a more in-depth explanation of [transformations](#transformations), [actions](#actions)
123 and [parallelism](#parallel-execution).
124 
125 ## <a name="crash-course"></a> Crash course
126 ### Filling a histogram
127 Let's start with a very common task: filling a histogram
128 ~~~{.cpp}
129 // Fill a TH1F with the "MET" branch
130 ROOT::Experimental::TDataFrame d("myTree", filePtr); // build a TDataFrame like you would build a TTreeReader
131 auto h = d.Histo1D("MET");
132 h->Draw();
133 ~~~
134 The first line creates a `TDataFrame` associated to the `TTree` "myTree". This tree has a branch named "MET".
135 
136 `Histo1D` is an action; it returns a smart pointer (a `TResultProxy` to be precise) to a `TH1F` histogram filled
137 with the `MET` of all events.
138 If the quantity stored in the branch is a collection, the histogram is filled with its elements.
139 
140 There are many other possible [actions](#overview), and all their results are wrapped in smart pointers; we'll see why
141 in a minute.
142 
143 ### Applying a filter
144 Let's now pretend we want to cut over the value of branch "MET" and count how many events pass this cut:
145 ~~~{.cpp}
146 // Select events with "MET" greater than 4., count events that passed the selection
147 auto metCut = [](double x) { return x > 4.; }; // a c++11 lambda function checking "x > 4"
148 ROOT::Experimental::TDataFrame d("myTree", filePtr);
149 auto c = d.Filter(metCut, {"MET"}).Count();
150 std::cout << *c << std::endl;
151 ~~~
152 `Filter` takes a function (a lambda in this example, but it can be any kind of function or even a functor class) and a
153 list of branch names. The filter function is applied to the specified branches for each event; it is required to return
154 a `bool` which signals whether the event passes the filter (`true`) or not (`false`). You can think of your data as
155 "flowing" through the chain of calls, being transformed, filtered and finally used to perform actions. Multiple `Filter`
156 calls can be chained one after another.
157 It is possible to specify filters as strings too. This snippet is analogous to the one above:
158 ~~~{.cpp}
159 ROOT::Experimental::TDataFrame d("myTree", filePtr);
160 auto c = d.Filter("MET > 4.").Count();
161 std::cout << *c << std::endl;
162 ~~~
163 Here the names of the branches used in the expression and their types are inferred automatically. The string must be
164 standard C++ and is just-in-time compiled by the ROOT interpreter, Cling.
165 
166 ### Running on a range of entries
167 It is sometimes necessary to limit the processing of the dataset to a range of entries. For this reason, the TDataFrame
168 offers the concept of ranges as a node of the TDataFrame graph: this means that filters, columns and actions can be hung
169 to it. If a range is specified after a filter, the range will act exclusively on the entries surviving the filter.
170 Here you can find some code using ranges:
171 ~~~{.cpp}
172 ROOT::Experimental::TDataFrame d("myTree", filePtr);
173 // This is how you can express a range of the first 30 entries
174 auto d_0_30 = d.Range(0, 30);
175 // This is how you pick all entries from 15 onwards
176 auto d_15_end = d.Range(15, 0);
177 // We can use a stride too, in this case we pick an event every 3
178 auto d_15_end_3 = d.Range(15, 0, 3);
179 ~~~
180 Ranges are not available when multi-threading is enabled.
181 
182 ### Creating a temporary column
183 Let's now consider the case in which "myTree" contains two quantities "x" and "y", but our analysis relies on a derived
184 quantity `z = sqrt(x*x + y*y)`.
185 Using the `Define` transformation, we can create a new column in the data-set containing the variable "z":
186 ~~~{.cpp}
187 auto sqrtSum = [](double x, double y) { return sqrt(x*x + y*y); };
188 auto zCut = [](double z) { return z > 0.; }
189 
190 ROOT::Experimental::TDataFrame d(treeName, filePtr);
191 auto zMean = d.Define("z", sqrtSum, {"x","y"})
192  .Filter(zCut, {"z"})
193  .Mean("z");
194 std::cout << *zMean << std::endl;
195 ~~~
196 `Define` creates the variable "z" by applying `sqrtSum` to "x" and "y". Later in the chain of calls we refer to
197 variables created with `Define` as if they were actual tree branches, but they are evaluated on the fly, once per
198 event. As with filters, `Define` calls can be chained with other transformations to create multiple temporary
199 columns.
200 As with filters, it is possible to specify new columns as strings too. This snippet is analogous to the one above:
201 ~~~{.cpp}
202 ROOT::Experimental::TDataFrame d(treeName, filePtr);
203 auto zMean = d.Define("z", "sqrt(x*x + y*y)")
204  .Filter("z > 0.")
205  .Mean("z");
206 std::cout << *zMean << std::endl;
207 ~~~
208 Again the names of the branches used in the expression and their types are inferred automatically. The string must be
209 standard C++ and is just-in-time compiled by the ROOT interpreter, Cling.
210 
211 ### Executing multiple actions
212 As a final example let us apply two different cuts on branch "MET" and fill two different histograms with the "pt\_v" of
213 the filtered events.
214 You should be able to easily understand what's happening:
215 ~~~{.cpp}
216 // fill two histograms with the results of two opposite cuts
217 auto isBig = [](double x) { return x > 10.; };
218 ROOT::Experimental::TDataFrame d(treeName, filePtr);
219 auto h1 = d.Filter(isBig, {"MET"}).Histo1D("pt_v");
220 auto h2 = d.Histo1D("pt_v");
221 h1->Draw(); // event loop is run once here
222 h2->Draw("SAME"); // no need to run the event loop again
223 ~~~
224 `TDataFrame` executes all above actions by **running the event-loop only once**. The trick is that actions are not
225 executed at the moment they are called, but they are **lazy**, i.e. delayed until the moment one of their results is
226 accessed through the smart pointer. At that time, the even loop is triggered and *all* results are produced
227 simultaneously.
228 
229 It is therefore good practice to declare all your filters and actions *before* accessing their results, allowing
230 `TDataFrame` to loop once and produce all results in one go.
231 
232 ### Going parallel
233 Let's say we would like to run the previous examples in parallel on several cores, dividing events fairly between cores.
234 The only modification required to the snippets would be the addition of this line *before* constructing the main
235 data-frame object:
236 ~~~{.cpp}
237 ROOT::EnableImplicitMT();
238 ~~~
239 Simple as that, enjoy your speed-up.
240 
241 ## <a name="more-features"></a>More features
242 Here is a list of the most important features that have been omitted in the "Crash course" for brevity's sake.
243 You don't need to read all these to start using `TDataFrame`, but they are useful to save typing time and runtime.
244 
245 ### Default branch lists
246 When constructing a `TDataFrame` object, it is possible to specify a **default branch list** for your analysis, in the
247 usual form of a list of strings representing branch names. The default branch list will be used as fallback whenever one
248 specific to the transformation/action is not present.
249 ~~~{.cpp}
250 // use "b1" and "b2" as default branches for `Filter`, `Define` and actions
251 ROOT::Experimental::TDataFrame d1(treeName, &file, {"b1","b2"});
252 // filter acts on default branch list, no need to specify it
253 auto h = d1.Filter([](int b1, int b2) { return b1 > b2; }).Histo1D("otherVar");
254 
255 // just one default branch this time
256 ROOT::Experimental::TDataFrame d2(treeName, &file, {"b1"});
257 // we can still specify non-default branch lists
258 // `Min` here can fall back to the default "b1"
259 auto min = d2.Filter([](double b2) { return b2 > 0; }, {"b2"}).Min();
260 ~~~
261 
262 ### Branch type guessing and explicit declaration of branch types
263 C++ is a statically typed language: all types must be known at compile-time. This includes the types of the `TTree`
264 branches we want to work on. For filters, temporary columns and some of the actions, **branch types are deduced from the
265 signature** of the relevant filter function/temporary column expression/action function:
266 ~~~{.cpp}
267 // here b1 is deduced to be `int` and b2 to be `double`
268 dataFrame.Filter([](int x, double y) { return x > 0 && y < 0.; }, {"b1", "b2"});
269 ~~~
270 If we specify an incorrect type for one of the branches, an exception with an informative message will be thrown at
271 runtime, when the branch value is actually read from the `TTree`: the implementation of `TDataFrame` allows the
272 detection of type mismatches. The same would happen if we swapped the order of "b1" and "b2" in the branch list passed
273 to `Filter`.
274 
275 Certain actions, on the other hand, do not take a function as argument (e.g. `Histo1D`), so we cannot deduce the type of
276 the branch at compile-time. In this case **`TDataFrame` tries to guess the type of the branch**, trying out the most
277 common ones and `std::vector` thereof. This is why we never needed to specify the branch types for all actions in the
278 above snippets.
279 
280 When the branch type is not a common one such as `int`, `double`, `char` or `float` it is therefore good practice to
281 specify it as a template parameter to the action itself, like this:
282 ~~~{.cpp}
283 dataFrame.Histo1D("b1"); // OK if b1 is a "common" type
284 dataFrame.Histo1D<Object_t>("myObject"); // OK, "myObject" is deduced to be of type `Object_t`
285 // dataFrame.Histo1D("myObject"); // THROWS an exception
286 ~~~
287 
288 ### Generic actions
289 `TDataFrame` strives to offer a comprehensive set of standard actions that can be performed on each event. At the same
290 time, it **allows users to execute arbitrary code (i.e. a generic action) inside the event loop** through the `Foreach`
291 and `ForeachSlot` actions.
292 
293 `Foreach(f, branchList)` takes a function `f` (lambda expression, free function, functor...) and a list of branches, and
294 executes `f` on those branches for each event. The function passed must return nothing (i.e. `void`). It can be used to
295 perform actions that are not already available in the interface. For example, the following snippet evaluates the root
296 mean square of branch "b":
297 ~~~{.cpp}
298 // Single-thread evaluation of RMS of branch "b" using Foreach
299 double sumSq = 0.;
300 unsigned int n = 0;
301 ROOT::Experimental::TDataFrame d("bTree", bFilePtr);
302 d.Foreach([&sumSq, &n](double b) { ++n; sumSq += b*b; }, {"b"});
303 std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
304 ~~~
305 When executing on multiple threads, users are responsible for the thread-safety of the expression passed to `Foreach`.
306 The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of `rms`; but
307 this is probably too much head-scratch for such a simple operation.
308 
309 `ForeachSlot` can help in this situation. It is an alternative version of `Foreach` for which the function takes an
310 additional parameter besides the branches it should be applied to: an `unsigned int slot` parameter, where `slot` is a
311 number indicating which thread (0, 1, 2 , ..., poolSize - 1) the function is being run in. We can take advantage of
312 `ForeachSlot` to evaluate a thread-safe root mean square of branch "b":
313 ~~~{.cpp}
314 // Thread-safe evaluation of RMS of branch "b" using ForeachSlot
315 ROOT::EnableImplicitMT();
316 unsigned int nSlots = ROOT::GetImplicitMTPoolSize();
317 std::vector<double> sumSqs(nSlots, 0.);
318 std::vector<unsigned int> ns(nSlots, 0);
319 
320 ROOT::Experimental::TDataFrame d("bTree", bFilePtr);
321 d.ForeachSlot([&sumSqs, &ns](unsigned int slot, double b) { sumSqs[slot] += b*b; ns[slot] += 1; }, {"b"});
322 double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
323 unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
324 std::cout << "rms of b: " << std::sqrt(sumSq / n) << std::endl;
325 ~~~
326 You see how we created one `double` variable for each thread in the pool, and later merged their results via
327 `std::accumulate`.
328 
329 ### Call graphs (storing and reusing sets of transformations)
330 **Sets of transformations can be stored as variables** and reused multiple times to create **call graphs** in which
331 several paths of filtering/creation of branches are executed simultaneously; we often refer to this as "storing the
332 state of the chain".
333 
334 This feature can be used, for example, to create a temporary column once and use it in several subsequent filters or
335 actions, or to apply a strict filter to the data-set *before* executing several other transformations and actions,
336 effectively reducing the amount of events processed.
337 
338 Let's try to make this clearer with a commented example:
339 ~~~{.cpp}
340 // build the data-frame and specify a default branch list
341 ROOT::Experimental::TDataFrame d(treeName, filePtr, {"var1", "var2", "var3"});
342 
343 // apply a cut and save the state of the chain
344 auto filtered = d.Filter(myBigCut);
345 
346 // plot branch "var1" at this point of the chain
347 auto h1 = filtered.Histo1D("var1");
348 
349 // create a new branch "vec" with a vector extracted from a complex object (only for filtered entries)
350 // and save the state of the chain
351 auto newBranchFiltered = filtered.Define("vec", [](const Obj& o) { return o.getVector(); }, {"obj"});
352 
353 // apply a cut and fill a histogram with "vec"
354 auto h2 = newBranchFiltered.Filter(cut1).Histo1D("vec");
355 
356 // apply a different cut and fill a new histogram
357 auto h3 = newBranchFiltered.Filter(cut2).Histo1D("vec");
358 
359 // Inspect results
360 h2->Draw(); // first access to an action result: run event-loop!
361 h3->Draw("SAME"); // event loop does not need to be run again here..
362 std::cout << "Entries in h1: " << h1->GetEntries() << std::endl; // ..or here
363 ~~~
364 `TDataFrame` detects when several actions use the same filter or the same temporary column, and **only evaluates each
365 filter or temporary column once per event**, regardless of how many times that result is used down the call graph.
366 Objects read from each branch are **built once and never copied**, for maximum efficiency.
367 When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated,
368 so it might be advisable to put the strictest filters first in the chain.
369 
370 ## <a name="transformations"></a>Transformations
371 ### Filters
372 A filter is defined through a call to `Filter(f, branchList)`. `f` can be a function, a lambda expression, a functor
373 class, or any other callable object. It must return a `bool` signalling whether the event has passed the selection
374 (`true`) or not (`false`). It must perform "read-only" actions on the branches, and should not have side-effects (e.g.
375 modification of an external or static variable) to ensure correct results when implicit multi-threading is active.
376 
377 `TDataFrame` only evaluates filters when necessary: if multiple filters are chained one after another, they are executed
378 in order and the first one returning `false` causes the event to be discarded and triggers the processing of the next
379 entry. If multiple actions or transformations depend on the same filter, that filter is not executed multiple times for
380 each entry: after the first access it simply serves a cached result.
381 
382 #### <a name="named-filters-and-cutflow-reports"></a>Named filters and cutflow reports
383 An optional string parameter `name` can be passed to the `Filter` method to create a **named filter**. Named filters
384 work as usual, but also keep track of how many entries they accept and reject.
385 
386 Statistics are retrieved through a call to the `Report` method:
387 
388 - when `Report` is called on the main `TDataFrame` object, it prints stats for all named filters declared up to that
389 point
390 - when called on a stored chain state (i.e. a chain/graph node), it prints stats for all named filters in the section of
391 the chain between the main `TDataFrame` and that node (included).
392 
393 Stats are printed in the same order as named filters have been added to the graph, and *refer to the latest event-loop*
394 that has been run using the relevant `TDataFrame`. If `Report` is called before the event-loop has been run at least
395 once, a run is triggered.
396 
397 ### Ranges
398 When `TDataFrame` is not being used in a multi-thread environment (i.e. no call to `EnableImplicitMT` was made),
399 `Range` transformations are available. These act very much like filters but instead of basing their decision on
400 a filter expression, they rely on `start`,`stop` and `stride` parameters.
401 
402 - `start`: number of entries that will be skipped before starting processing again
403 - `stop`: maximum number of entries that will be processed
404 - `stride`: only process one entry every `stride` entries
405 
406 The actual number of entries processed downstream of a `Range` node will be `(stop - start)/stride` (or less if less
407 entries than that are available).
408 
409 Note that ranges act "locally", not based on the global entry count: `Range(10,50)` means "skip the first 10 entries
410 *that reach this node*, let the next 40 entries pass, then stop processing". If a range node hangs from a filter node,
411 and the range has a `start` parameter of 10, that means the range will skip the first 10 entries *that pass the
412 preceding filter*.
413 
414 Ranges allow "early quitting": if all branches of execution of a functional graph reached their `stop` value of
415 processed entries, the event-loop is immediately interrupted. This is useful for debugging and initial explorations.
416 
417 ### Temporary columns
418 Temporary columns are created by invoking `Define(name, f, branchList)`. As usual, `f` can be any callable object
419 (function, lambda expression, functor class...); it takes the values of the branches listed in `branchList` (a list of
420 strings) as parameters, in the same order as they are listed in `branchList`. `f` must return the value that will be
421 assigned to the temporary column.
422 
423 A new variable is created called `name`, accessible as if it was contained in the dataset from subsequent
424 transformations/actions.
425 
426 Use cases include:
427 - caching the results of complex calculations for easy and efficient multiple access
428 - extraction of quantities of interest from complex objects
429 - branch aliasing, i.e. changing the name of a branch
430 
431 An exception is thrown if the `name` of the new branch is already in use for another branch in the `TTree`.
432 
433 It is also possible to specify the quantity to be stored in the new temporary column as a C++ expression with the method
434 `Define(name, expression)`. For example this invocation
435 
436 ~~~{.cpp}
437 tdf.Define("pt", "sqrt(px*px + py*py)");
438 ~~~
439 
440 will create a new column called "pt" the value of which is calculated starting from the branches px and py. The system
441 builds a just-in-time compiled function starting from the expression after having deduced the list of necessary branches
442 from the names of the variables specified by the user.
443 
444 ## <a name="actions"></a>Actions
445 ### Instant and lazy actions
446 Actions can be **instant** or **lazy**. Instant actions are executed as soon as they are called, while lazy actions are
447 executed whenever the object they return is accessed for the first time. As a rule of thumb, actions with a return value
448 are lazy, the others are instant.
449 
450 ### Overview
451 Here is a quick overview of what actions are present and what they do. Each one is described in more detail in the
452 reference guide.
453 
454 In the following, whenever we say an action "returns" something, we always mean it returns a smart pointer to it. Also
455 note that all actions are only executed for events that pass all preceding filters.
456 
457 | **Lazy actions** | **Description** |
458 |------------------|-----------------|
459 | Count | Return the number of events processed. |
460 | Fill | Fill a user-defined object with the values of the specified branches, as if by calling `Obj.Fill(branch1, branch2, ...). |
461 | Histo{1D,2D,3D} | Fill a {one,two,three}-dimensional histogram with the processed branch values. |
462 | Max | Return the maximum of processed branch values. |
463 | Mean | Return the mean of processed branch values. |
464 | Min | Return the minimum of processed branch values. |
465 | Profile{1D,2D} | Fill a {one,two}-dimensional profile with the branch values that passed all filters. |
466 | 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 branch. Return the final result of the reduction operation. An optional parameter allows initialization of the result object to non-default values. |
467 | Take | Build a collection of values of a branch. |
468 
469 | **Instant actions** | **Description** |
470 |---------------------|-----------------|
471 | 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. |
472 | 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 `TDataFrame` after `ROOT::EnableImplicitMT()`. `ForeachSlot` works just as well with single-thread execution: in that case `slot` will always be `0`. |
473 | Snapshot | Writes on disk a dataset made of the selected columns and entries passing the filters (if any). |
474 
475 | **Queries** | **Description** |
476 |-----------|-----------------|
477 | Report | This is not properly an action, since when `Report` is called it does not book an operation to be performed on each entry. Instead, it interrogates the data-frame directly to print a cutflow report, i.e. 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. |
478 
479 ## <a name="parallel-execution"></a>Parallel execution
480 As pointed out before in this document, `TDataFrame` can transparently perform multi-threaded event loops to speed up
481 the execution of its actions. Users only have to call `ROOT::EnableImplicitMT()` *before* constructing the `TDataFrame`
482 object to indicate that it should take advantage of a pool of worker threads. **Each worker thread processes a distinct
483 subset of entries**, and their partial results are merged before returning the final values to the user.
484 
485 ### Thread safety
486 `Filter` and `Define` transformations should be inherently thread-safe: they have no side-effects and are not
487 dependent on global state.
488 Most `Filter`/`Define` functions will in fact be pure in the functional programming sense.
489 All actions are built to be thread-safe with the exception of `Foreach`, in which case users are responsible of
490 thread-safety, see [here](#generic-actions).
491 
492 <a name="reference"></a>
493 */
494 
495 ////////////////////////////////////////////////////////////////////////////
496 /// \brief Build the dataframe
497 /// \param[in] treeName Name of the tree contained in the directory
498 /// \param[in] dirPtr TDirectory where the tree is stored, e.g. a TFile.
499 /// \param[in] defaultBranches Collection of default branches.
500 ///
501 /// The default branches are looked at in case no branch is specified in the
502 /// booking of actions or transformations.
503 /// See TInterface for the documentation of the
504 /// methods available.
505 TDataFrame::TDataFrame(std::string_view treeName, TDirectory *dirPtr, const ColumnNames_t &defaultBranches)
506  : TInterface<TDFDetail::TLoopManager>(std::make_shared<TDFDetail::TLoopManager>(nullptr, defaultBranches))
507 {
508  const std::string treeNameInt(treeName);
509  if (!dirPtr) {
510  auto msg = "Invalid TDirectory!";
511  throw std::runtime_error(msg);
512  }
513  auto tree = static_cast<TTree *>(dirPtr->Get(treeNameInt.c_str()));
514  if (!tree) {
515  auto msg = "Tree \"" + treeNameInt + "\" cannot be found!";
516  throw std::runtime_error(msg);
517  }
518  fProxiedPtr->SetTree(std::shared_ptr<TTree>(tree, [](TTree *) {}));
519 }
520 
521 ////////////////////////////////////////////////////////////////////////////
522 /// \brief Build the dataframe
523 /// \param[in] treeName Name of the tree contained in the directory
524 /// \param[in] filenameglob TDirectory where the tree is stored, e.g. a TFile.
525 /// \param[in] defaultBranches Collection of default branches.
526 ///
527 /// The default branches are looked at in case no branch is specified in the
528 /// booking of actions or transformations.
529 /// See TInterface for the documentation of the
530 /// methods available.
531 TDataFrame::TDataFrame(std::string_view treeName, std::string_view filenameglob,
532  const ColumnNames_t &defaultBranches)
533  : TInterface<TDFDetail::TLoopManager>(std::make_shared<TDFDetail::TLoopManager>(nullptr, defaultBranches))
534 {
535  const std::string treeNameInt(treeName);
536  const std::string filenameglobInt(filenameglob);
537  auto chain = new TChain(treeNameInt.c_str());
538  chain->Add(filenameglobInt.c_str());
539  fProxiedPtr->SetTree(std::shared_ptr<TTree>(static_cast<TTree *>(chain)));
540 }
541 
542 ////////////////////////////////////////////////////////////////////////////
543 /// \brief Build the dataframe
544 /// \param[in] tree The tree or chain to be studied.
545 /// \param[in] defaultBranches Collection of default branches.
546 ///
547 /// The default branches are looked at in case no branch is specified in the
548 /// booking of actions or transformations.
549 /// See TInterface for the documentation of the
550 /// methods available.
551 TDataFrame::TDataFrame(TTree &tree, const ColumnNames_t &defaultBranches)
552  : TInterface<TDFDetail::TLoopManager>(std::make_shared<TDFDetail::TLoopManager>(&tree, defaultBranches))
553 {
554 }
555 
556 //////////////////////////////////////////////////////////////////////////
557 /// \brief Build the dataframe
558 /// \param[in] numEntries The number of entries to generate.
559 ///
560 /// An empty-source dataframe constructed with a number of entries will
561 /// generate those entries on the fly when some action is triggered,
562 /// and it will do so for all the previously-defined temporary branches.
564  : TInterface<TDFDetail::TLoopManager>(std::make_shared<TDFDetail::TLoopManager>(numEntries))
565 {
566 }
TDataFrame(std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultBranches={})
Build the dataframe.
Definition: TDataFrame.cxx:531
long long Long64_t
Definition: RtypesCore.h:69
TDFDetail::ColumnNames_t ColumnNames_t
Definition: TDataFrame.hxx:37
STL namespace.
std::shared_ptr< ToContentType_t< T > > Get(const std::string &name)
Get the object for a key.
Definition: TDirectory.hxx:149
std::shared_ptr< TDFDetail::TLoopManager > fProxiedPtr
Key/value store of objects.
Definition: TDirectory.hxx:68
A chain is a collection of files containing TTree objects.
Definition: TChain.h:33
Definition: tree.py:1
virtual Int_t Add(TChain *chain)
Add all files referenced by the passed chain to this chain.
Definition: TChain.cxx:220