Logo ROOT  
Reference Guide
df001_introduction.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -nodraw
4/// This tutorial illustrates the basic features of the RDataFrame class,
5/// a utility which allows to interact with data stored in TTrees following
6/// a functional-chain like approach.
7///
8/// \macro_code
9/// \macro_output
10///
11/// \date December 2016
12/// \author Enrico Guiraud
13
14// ## Preparation
15
16// A simple helper function to fill a test tree: this makes the example
17// stand-alone.
18void fill_tree(const char *treeName, const char *fileName)
19{
21 int i(0);
22 d.Define("b1", [&i]() { return (double)i; })
23 .Define("b2",
24 [&i]() {
25 auto j = i * i;
26 ++i;
27 return j;
28 })
29 .Snapshot(treeName, fileName);
30}
31
33{
34
35 // We prepare an input tree to run on
36 auto fileName = "df001_introduction.root";
37 auto treeName = "myTree";
38 fill_tree(treeName, fileName);
39
40 // We read the tree from the file and create a RDataFrame, a class that
41 // allows us to interact with the data contained in the tree.
42 // We select a default column, a *branch* to adopt ROOT jargon, which will
43 // be looked at if none is specified by the user when dealing with filters
44 // and actions.
45 ROOT::RDataFrame d(treeName, fileName, {"b1"});
46
47 // ## Operations on the dataframe
48 // We now review some *actions* which can be performed on the data frame.
49 // All actions but ForEach return a TActionResultPtr<T>. The series of
50 // operations on the data frame is not executed until one of those pointers
51 // is accessed. If the Foreach action is invoked, the execution is immediate.
52 // But first of all, let us we define now our cut-flow with two lambda
53 // functions. We can use free functions too.
54 auto cutb1 = [](double b1) { return b1 < 5.; };
55 auto cutb1b2 = [](int b2, double b1) { return b2 % 2 && b1 < 4.; };
56
57 // ### `Count` action
58 // The `Count` allows to retrieve the number of the entries that passed the
59 // filters. Here we show how the automatic selection of the column kicks
60 // in in case the user specifies none.
61 auto entries1 = d.Filter(cutb1) // <- no column name specified here!
62 .Filter(cutb1b2, {"b2", "b1"})
63 .Count();
64
65 std::cout << *entries1 << " entries passed all filters" << std::endl;
66
67 // Filters can be expressed as strings. The content must be C++ code. The
68 // name of the variables must be the name of the branches. The code is
69 // just in time compiled.
70 auto entries2 = d.Filter("b1 < 5.").Count();
71 std::cout << *entries2 << " entries passed the string filter" << std::endl;
72
73 // ### `Min`, `Max` and `Mean` actions
74 // These actions allow to retrieve statistical information about the entries
75 // passing the cuts, if any.
76 auto b1b2_cut = d.Filter(cutb1b2, {"b2", "b1"});
77 auto minVal = b1b2_cut.Min();
78 auto maxVal = b1b2_cut.Max();
79 auto meanVal = b1b2_cut.Mean();
80 auto nonDefmeanVal = b1b2_cut.Mean("b2"); // <- Column is not the default
81 std::cout << "The mean is always included between the min and the max: " << *minVal << " <= " << *meanVal
82 << " <= " << *maxVal << std::endl;
83
84 // ### `Take` action
85 // The `Take` action allows to retrieve all values of the variable stored in a
86 // particular column that passed filters we specified. The values are stored
87 // in a list by default, but other collections can be chosen.
88 auto b1_cut = d.Filter(cutb1);
89 auto b1Vec = b1_cut.Take<double>();
90 auto b1List = b1_cut.Take<double, std::list<double>>();
91
92 std::cout << "Selected b1 entries" << std::endl;
93 for (auto b1_entry : *b1List)
94 std::cout << b1_entry << " ";
95 std::cout << std::endl;
96 auto b1VecCl = ROOT::GetClass(b1Vec.GetPtr());
97 std::cout << "The type of b1Vec is " << b1VecCl->GetName() << std::endl;
98
99 // ### `Histo1D` action
100 // The `Histo1D` action allows to fill an histogram. It returns a TH1D filled
101 // with values of the column that passed the filters. For the most common
102 // types, the type of the values stored in the column is automatically
103 // guessed.
104 auto hist = d.Filter(cutb1).Histo1D();
105 std::cout << "Filled h " << hist->GetEntries() << " times, mean: " << hist->GetMean() << std::endl;
106
107 // ### `Foreach` action
108 // The most generic action of all: an operation is applied to all entries.
109 // In this case we fill a histogram. In some sense this is a violation of a
110 // purely functional paradigm - C++ allows to do that.
111 TH1F h("h", "h", 12, -1, 11);
112 d.Filter([](int b2) { return b2 % 2 == 0; }, {"b2"}).Foreach([&h](double b1) { h.Fill(b1); });
113
114 std::cout << "Filled h with " << h.GetEntries() << " entries" << std::endl;
115
116 // ## Express your chain of operations with clarity!
117 // We are discussing an example here but it is not hard to imagine much more
118 // complex pipelines of actions acting on data. Those might require code
119 // which is well organised, for example allowing to conditionally add filters
120 // or again to clearly separate filters and actions without the need of
121 // writing the entire pipeline on one line. This can be easily achieved.
122 // We'll show this re-working the `Count` example:
123 auto cutb1_result = d.Filter(cutb1);
124 auto cutb1b2_result = d.Filter(cutb1b2, {"b2", "b1"});
125 auto cutb1_cutb1b2_result = cutb1_result.Filter(cutb1b2, {"b2", "b1"});
126 // Now we want to count:
127 auto evts_cutb1_result = cutb1_result.Count();
128 auto evts_cutb1b2_result = cutb1b2_result.Count();
129 auto evts_cutb1_cutb1b2_result = cutb1_cutb1b2_result.Count();
130
131 std::cout << "Events passing cutb1: " << *evts_cutb1_result << std::endl
132 << "Events passing cutb1b2: " << *evts_cutb1b2_result << std::endl
133 << "Events passing both: " << *evts_cutb1_cutb1b2_result << std::endl;
134
135 // ## Calculating quantities starting from existing columns
136 // Often, operations need to be carried out on quantities calculated starting
137 // from the ones present in the columns. We'll create in this example a third
138 // column the values of which are the sum of the *b1* and *b2* ones, entry by
139 // entry. The way in which the new quantity is defined is via a runable.
140 // It is important to note two aspects at this point:
141 // - The value is created on the fly only if the entry passed the existing
142 // filters.
143 // - The newly created column behaves as the one present on the file on disk.
144 // - The operation creates a new value, without modifying anything. De facto,
145 // this is like having a general container at disposal able to accommodate
146 // any value of any type.
147 // Let's dive in an example:
148 auto entries_sum = d.Define("sum", [](double b1, int b2) { return b2 + b1; }, {"b1", "b2"})
149 .Filter([](double sum) { return sum > 4.2; }, {"sum"})
150 .Count();
151 std::cout << *entries_sum << std::endl;
152
153 // Additional columns can be expressed as strings. The content must be C++
154 // code. The name of the variables must be the name of the branches. The code
155 // is just in time compiled.
156 auto entries_sum2 = d.Define("sum2", "b1 + b2").Filter("sum2 > 4.2").Count();
157 std::cout << *entries_sum2 << std::endl;
158
159 // It is possible at any moment to read the entry number and the processing
160 // slot number. The latter may change when implicit multithreading is active.
161 // The special columns which provide the entry number and the slot index are
162 // called "rdfentry_" and "rdfslot_" respectively. Their types are an unsigned
163 // 64 bit integer and an unsigned integer.
164 auto printEntrySlot = [](ULong64_t iEntry, unsigned int slot) {
165 std::cout << "Entry: " << iEntry << " Slot: " << slot << std::endl;
166 };
167 d.Foreach(printEntrySlot, {"rdfentry_", "rdfslot_"});
168
169 return 0;
170}
double
Definition: Converters.cxx:921
#define d(i)
Definition: RSha256.hxx:102
#define h(i)
Definition: RSha256.hxx:106
unsigned long long ULong64_t
Definition: RtypesCore.h:72
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition: RVec.hxx:939
TClass * GetClass(T *)
Definition: TClass.h:658
static long int sum(long int i)
Definition: Factory.cxx:2275