25#include <nlohmann/json.hpp>
1344 void Exec(unsigned int slot)
1346 fPerThreadResults[slot]++;
1349 // Called at the end of the event loop.
1352 *fFinalResult = std::accumulate(fPerThreadResults.begin(), fPerThreadResults.end(), 0);
1355 // Called by RDataFrame to retrieve the name of this action.
1356 std::string GetActionName() const { return "MyCounter"; }
1360 ROOT::RDataFrame df(10);
1361 ROOT::RDF::RResultPtr<int> resultPtr = df.Book<>(MyCounter{df.GetNSlots()}, {});
1362 // The GetValue call triggers the event loop
1363 std::cout << "Number of processed entries: " << resultPtr.GetValue() << std::endl;
1367See the Book() method for more information and [this tutorial](https://root.cern/doc/master/df018__customActions_8C.html)
1368for a more complete example.
1370#### Injecting arbitrary code in the event loop with Foreach() and ForeachSlot()
1372Foreach() takes a callable (lambda expression, free function, functor...) and a list of columns and
1373executes the callable on the values of those columns for each event that passes all upstream selections.
1374It can be used to perform actions that are not already available in the interface. For example, the following snippet
1375evaluates the root mean square of column "x":
1377// Single-thread evaluation of RMS of column "x" using Foreach
1380df.Foreach([&sumSq, &n](double x) { ++n; sumSq += x*x; }, {"x"});
1381std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1383In multi-thread runs, users are responsible for the thread-safety of the expression passed to Foreach():
1384thread will execute the expression concurrently.
1385The code above would need to employ some resource protection mechanism to ensure non-concurrent writing of `rms`; but
1386this is probably too much head-scratch for such a simple operation.
1388ForeachSlot() can help in this situation. It is an alternative version of Foreach() for which the function takes an
1389additional "processing slot" parameter besides the columns it should be applied to. RDataFrame
1390guarantees that ForeachSlot() will invoke the user expression with different `slot` parameters for different concurrent
1391executions (see [Special helper columns: rdfentry_ and rdfslot_](\ref helper-cols) for more information on the slot parameter).
1392We can take advantage of ForeachSlot() to evaluate a thread-safe root mean square of column "x":
1394// Thread-safe evaluation of RMS of column "x" using ForeachSlot
1395ROOT::EnableImplicitMT();
1396const unsigned int nSlots = df.GetNSlots();
1397std::vector<double> sumSqs(nSlots, 0.);
1398std::vector<unsigned int> ns(nSlots, 0);
1400df.ForeachSlot([&sumSqs, &ns](unsigned int slot, double x) { sumSqs[slot] += x*x; ns[slot] += 1; }, {"x"});
1401double sumSq = std::accumulate(sumSqs.begin(), sumSqs.end(), 0.); // sum all squares
1402unsigned int n = std::accumulate(ns.begin(), ns.end(), 0); // sum all counts
1403std::cout << "rms of x: " << std::sqrt(sumSq / n) << std::endl;
1405Notice how we created one `double` variable for each processing slot and later merged their results via `std::accumulate`.
1409### Dataset joins with friend trees
1411Vertically concatenating multiple trees that have the same columns (creating a logical dataset with the same columns and
1412more rows) is trivial in RDataFrame: just pass the tree name and a list of file names to RDataFrame's constructor, or create a TChain
1413out of the desired trees and pass that to RDataFrame.
1415Horizontal concatenations of trees or chains (creating a logical dataset with the same number of rows and the union of the
1416columns of multiple trees) leverages TTree's "friend" mechanism.
1418Simple joins of trees that do not have the same number of rows are also possible with indexed friend trees (see below).
1420To use friend trees in RDataFrame, set up trees with the appropriate relationships and then instantiate an RDataFrame
1426main.AddFriend(&friend, "myFriend");
1429auto df2 = df.Filter("myFriend.MyCol == 42");
1432The same applies for TChains. Columns coming from the friend trees can be referred to by their full name, like in the example above,
1433or the friend tree name can be omitted in case the column name is not ambiguous (e.g. "MyCol" could be used instead of
1434"myFriend.MyCol" in the example above if there is no column "MyCol" in the main tree).
1436\note A common source of confusion is that trees that are written out from a multi-thread Snapshot() call will have their
1437 entries (block-wise) shuffled with respect to the original tree. Such trees cannot be used as friends of the original
1438 one: rows will be mismatched.
1440Indexed friend trees provide a way to perform simple joins of multiple trees over a common column.
1441When a certain entry in the main tree (or chain) is loaded, the friend trees (or chains) will then load an entry where the
1442"index" columns have a value identical to the one in the main one. For example, in Python:
1448# If a friend tree has an index on `commonColumn`, when the main tree loads
1449# a given row, it also loads the row of the friend tree that has the same
1450# value of `commonColumn`
1451aux_tree.BuildIndex("commonColumn")
1453mainTree.AddFriend(aux_tree)
1455df = ROOT.RDataFrame(mainTree)
1458RDataFrame supports indexed friend TTrees from ROOT v6.24 in single-thread mode and from v6.28/02 in multi-thread mode.
1460\anchor other-file-formats
1461### Reading data formats other than ROOT trees
1462RDataFrame can be interfaced with RDataSources. The ROOT::RDF::RDataSource interface defines an API that RDataFrame can use to read arbitrary columnar data formats.
1464RDataFrame calls into concrete RDataSource implementations to retrieve information about the data, retrieve (thread-local) readers or "cursors" for selected columns
1465and to advance the readers to the desired data entry.
1466Some predefined RDataSources are natively provided by ROOT such as the ROOT::RDF::RCsvDS which allows to read comma separated files:
1468auto tdf = ROOT::RDF::FromCSV("MuRun2010B.csv");
1469auto filteredEvents =
1470 tdf.Filter("Q1 * Q2 == -1")
1471 .Define("m", "sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))");
1472auto h = filteredEvents.Histo1D("m");
1476See also FromNumpy (Python-only), FromRNTuple(), FromArrow(), FromSqlite().
1479### Computation graphs (storing and reusing sets of transformations)
1481As we saw, transformed dataframes can be stored as variables and reused multiple times to create modified versions of the dataset. This implicitly defines a **computation graph** in which
1482several paths of filtering/creation of columns are executed simultaneously, and finally aggregated results are produced.
1484RDataFrame detects when several actions use the same filter or the same defined column, and **only evaluates each
1485filter or defined column once per event**, regardless of how many times that result is used down the computation graph.
1486Objects read from each column are **built once and never copied**, for maximum efficiency.
1487When "upstream" filters are not passed, subsequent filters, temporary column expressions and actions are not evaluated,
1488so it might be advisable to put the strictest filters first in the graph.
1490\anchor representgraph
1491### Visualizing the computation graph
1492It is possible to print the computation graph from any node to obtain a [DOT (graphviz)](https://en.wikipedia.org/wiki/DOT_(graph_description_language)) representation either on the standard output
1495Invoking the function ROOT::RDF::SaveGraph() on any node that is not the head node, the computation graph of the branch
1496the node belongs to is printed. By using the head node, the entire computation graph is printed.
1498Following there is an example of usage:
1500// First, a sample computational graph is built
1501ROOT::RDataFrame df("tree", "f.root");
1503auto df2 = df.Define("x", []() { return 1; })
1504 .Filter("col0 % 1 == col0")
1505 .Filter([](int b1) { return b1 <2; }, {"cut1"})
1506 .Define("y", []() { return 1; });
1508auto count = df2.Count();
1510// Prints the graph to the rd1.dot file in the current directory
1511ROOT::RDF::SaveGraph(df, "./mydot.dot");
1512// Prints the graph to standard output
1513ROOT::RDF::SaveGraph(df);
1516The generated graph can be rendered using one of the graphviz filters, e.g. `dot`. For instance, the image below can be generated with the following command:
1518$ dot -Tpng computation_graph.dot -ocomputation_graph.png
1521\image html RDF_Graph2.png
1524### Activating RDataFrame execution logs
1526RDataFrame has experimental support for verbose logging of the event loop runtimes and other interesting related information. It is activated as follows:
1528#include <ROOT/RLogger.hxx>
1530// this increases RDF's verbosity level as long as the `verbosity` variable is in scope
1531auto verbosity = ROOT::Experimental::RLogScopedVerbosity(ROOT::Detail::RDF::RDFLogChannel(), ROOT::Experimental::ELogLevel::kInfo);
1538verbosity = ROOT.Experimental.RLogScopedVerbosity(ROOT.Detail.RDF.RDFLogChannel(), ROOT.Experimental.ELogLevel.kInfo)
1541More information (e.g. start and end of each multi-thread task) is printed using `ELogLevel.kDebug` and even more
1542(e.g. a full dump of the generated code that RDataFrame just-in-time-compiles) using `ELogLevel.kDebug+10`.
1544\anchor rdf-from-spec
1545### Creating an RDataFrame from a dataset specification file
1547RDataFrame can be created using a dataset specification JSON file:
1552df = ROOT.RDF.Experimental.FromSpec("spec.json")
1555The input dataset specification JSON file needs to be provided by the user and it describes all necessary samples and
1556their associated metadata information. The main required key is the "samples" (at least one sample is needed) and the
1557required sub-keys for each sample are "trees" and "files". Additionally, one can specify a metadata dictionary for each
1558sample in the "metadata" key.
1560A simple example for the formatting of the specification in the JSON file is the following:
1566 "trees": ["tree1", "tree2"],
1567 "files": ["file1.root", "file2.root"],
1571 "sample_category" = "data"
1575 "trees": ["tree3", "tree4"],
1576 "files": ["file3.root", "file4.root"],
1580 "sample_category" = "MC_background"
1587The metadata information from the specification file can be then accessed using the DefinePerSample function.
1588For example, to access luminosity information (stored as a double):
1591df.DefinePerSample("lumi", 'rdfsampleinfo_.GetD("lumi")')
1594or sample_category information (stored as a string):
1597df.DefinePerSample("sample_category", 'rdfsampleinfo_.GetS("sample_category")')
1600or directly the filename:
1603df.DefinePerSample("name", "rdfsampleinfo_.GetSampleName()")
1606An example implementation of the "FromSpec" method is available in tutorial: df106_HiggstoFourLeptons.py, which also
1607provides a corresponding exemplary JSON file for the dataset specification.
1610### Adding a progress bar
1612A progress bar showing the processed event statistics can be added to any RDataFrame program.
1613The event statistics include elapsed time, currently processed file, currently processed events, the rate of event processing
1614and an estimated remaining time (per file being processed). It is recorded and printed in the terminal every m events and every
1615n seconds (by default m = 1000 and n = 1). The ProgressBar can be also added when the multithread (MT) mode is enabled.
1617ProgressBar is added after creating the dataframe object (df):
1619ROOT::RDataFrame df("tree", "file.root");
1620ROOT::RDF::Experimental::AddProgressBar(df);
1623Alternatively, RDataFrame can be cast to an RNode first, giving the user more flexibility
1624For example, it can be called at any computational node, such as Filter or Define, not only the head node,
1625with no change to the ProgressBar function itself (please see the [Python interface](classROOT_1_1RDataFrame.html#python)
1626section for appropriate usage in Python):
1628ROOT::RDataFrame df("tree", "file.root");
1629auto df_1 = ROOT::RDF::RNode(df.Filter("x>1"));
1630ROOT::RDF::Experimental::AddProgressBar(df_1);
1632Examples of implemented progress bars can be seen by running [Higgs to Four Lepton tutorial](https://root.cern/doc/master/df106__HiggsToFourLeptons_8py_source.html) and [Dimuon tutorial](https://root.cern/doc/master/df102__NanoAODDimuonAnalysis_8C.html).
1634\anchor missing-values
1635### Working with missing values in the dataset
1637In certain situations a dataset might be missing one or more values at one or
1638more of its entries. For example:
1640- If the dataset is composed of multiple files and one or more files is
1641 missing one or more columns required by the analysis.
1642- When joining different datasets horizontally according to some index value
1643 (e.g. the event number), if the index does not find a match in one or more
1644 other datasets for a certain entry.
1646For example, suppose that column "y" does not have a value for entry 42:
1656If the RDataFrame application reads that column, for example if a Take() action
1657was requested, the default behaviour is to throw an exception indicating
1658that that column is missing an entry.
1660The following paragraphs discuss the functionalities provided by RDataFrame to
1661work with missing values in the dataset.
1663#### FilterAvailable and FilterMissing
1665FilterAvailable and FilterMissing are specialized RDataFrame Filter operations.
1666They take as input argument the name of a column of the dataset to watch for
1667missing values. Like Filter, they will either keep or discard an entire entry
1668based on whether a condition returns true or false. Specifically:
1670- FilterAvailable: the condition is whether the value of the column is present.
1671 If so, the entry is kept. Otherwise if the value is missing the entry is
1673- FilterMissing: the condition is whether the value of the column is missing. If
1674 so, the entry is kept. Otherwise if the value is present the entry is
1678df = ROOT.RDataFrame(dataset)
1680# Anytime an entry from "col" is missing, the entire entry will be filtered out
1681df_available = df.FilterAvailable("col")
1682df_available = df_available.Define("twice", "col * 2")
1684# Conversely, if we want to select the entries for which the column has missing
1685# values, we do the following
1686df_missingcol = df.FilterMissing("col")
1687# Following operations in the same branch of the computation graph clearly
1688# cannot access that same column, since there would be no value to read
1689df_missingcol = df_missingcol.Define("observable", "othercolumn * 2")
1693ROOT::RDataFrame df{dataset};
1695// Anytime an entry from "col" is missing, the entire entry will be filtered out
1696auto df_available = df.FilterAvailable("col");
1697auto df_twicecol = df_available.Define("twice", "col * 2");
1699// Conversely, if we want to select the entries for which the column has missing
1700// values, we do the following
1701auto df_missingcol = df.FilterMissing("col");
1702// Following operations in the same branch of the computation graph clearly
1703// cannot access that same column, since there would be no value to read
1704auto df_observable = df_missingcol.Define("observable", "othercolumn * 2");
1709DefaultValueFor creates a node of the computation graph which just forwards the
1710values of the columns necessary for other downstream nodes, when they are
1711available. In case a value of the input column passed to this function is not
1712available, the node will provide the default value passed to this function call
1716df = ROOT.RDataFrame(dataset)
1717# Anytime an entry from "col" is missing, the value will be the default one
1718default_value = ... # Some sensible default value here
1719df = df.DefaultValueFor("col", default_value)
1720df = df.Define("twice", "col * 2")
1724ROOT::RDataFrame df{dataset};
1725// Anytime an entry from "col" is missing, the value will be the default one
1726constexpr auto default_value = ... // Some sensible default value here
1727auto df_default = df.DefaultValueFor("col", default_value);
1728auto df_col = df_default.Define("twice", "col * 2");
1731#### Mixing different strategies to work with missing values in the same RDataFrame
1733All the operations presented above only act on the particular branch of the
1734computation graph where they are called, so that different results can be
1735obtained by mixing and matching the filtering or providing a default value
1739df = ROOT.RDataFrame(dataset)
1740# Anytime an entry from "col" is missing, the value will be the default one
1741default_value = ... # Some sensible default value here
1742df_default = df.DefaultValueFor("col", default_value).Define("twice", "col * 2")
1743df_filtered = df.FilterAvailable("col").Define("twice", "col * 2")
1745# Same number of total entries as the input dataset, with defaulted values
1746df_default.Display(["twice"]).Print()
1747# Only keep the entries where "col" has values
1748df_filtered.Display(["twice"]).Print()
1752ROOT::RDataFrame df{dataset};
1754// Anytime an entry from "col" is missing, the value will be the default one
1755constexpr auto default_value = ... // Some sensible default value here
1756auto df_default = df.DefaultValueFor("col", default_value).Define("twice", "col * 2");
1757auto df_filtered = df.FilterAvailable("col").Define("twice", "col * 2");
1759// Same number of total entries as the input dataset, with defaulted values
1760df_default.Display({"twice"})->Print();
1761// Only keep the entries where "col" has values
1762df_filtered.Display({"twice"})->Print();
1765#### Further considerations
1767Note that working with missing values is currently supported with a TTree-based
1768data source. Support of this functionality for other data sources may come in
1792 auto msg =
"Invalid TDirectory!";
1793 throw std::runtime_error(
msg);
1799 throw std::runtime_error(
msg);
1801 GetProxiedPtr()->SetTree(std::shared_ptr<TTree>(tree, [](
TTree *) {}));
1928namespace Experimental {
1962 const nlohmann::ordered_json
fullData = nlohmann::ordered_json::parse(std::ifstream(
jsonFile));
1964 throw std::runtime_error(
1965 R
"(The input specification does not contain any samples. Please provide the samples in the specification like:
1969 "trees": ["tree1", "tree2"],
1970 "files": ["file1.root", "file2.root"],
1971 "metadata": {"lumi": 1.0, }
1974 "trees": ["tree3", "tree4"],
1975 "files": ["file3.root", "file4.root"],
1976 "metadata": {"lumi": 0.5, }
1989 if (!
sample.contains(
"trees")) {
1990 throw std::runtime_error(
"A list of tree names must be provided for sample " +
sampleName +
".");
1993 if (!
sample.contains(
"files")) {
1994 throw std::runtime_error(
"A list of files must be provided for sample " +
sampleName +
".");
1997 if (!
sample.contains(
"metadata")) {
2001 for (
const auto &metadata :
sample[
"metadata"].items()) {
2002 const auto &val = metadata.value();
2003 if (val.is_string())
2004 m.Add(metadata.key(), val.get<std::string>());
2005 else if (val.is_number_integer())
2006 m.Add(metadata.key(), val.get<
int>());
2007 else if (val.is_number_float())
2008 m.Add(metadata.key(), val.get<
double>());
2010 throw std::logic_error(
"The metadata keys can only be of type [string|int|double].");
2015 if (
fullData.contains(
"friends")) {
2018 std::vector<std::string>
trees =
friends.value()[
"trees"];
2019 std::vector<std::string>
files =
friends.value()[
"files"];
2021 throw std::runtime_error(
"Mismatch between trees and files in a friend.");
2029 if (
range.size() == 1)
2031 else if (
range.size() == 2)
2056 throw std::runtime_error(
"Cannot print information about this RDataFrame, "
2057 "it was not properly created. It must be discarded.");
2060 auto defCols =
lm->GetDefaultColumnNames();
2062 std::ostringstream
ret;
2069 ret <<
"\nDefault columns:\n";
2071 ret <<
" - " << col <<
"\n";
2076 ret <<
"A data frame associated to the data source \"" << cling::printValue(
ds) <<
"\"";
2078 ret <<
"An empty data frame that will create " <<
lm->GetNEmptyEntries() <<
" entries\n";
unsigned long long ULong64_t
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
The head node of a RDF computation graph.
The dataset specification for RDataFrame.
Class representing a sample which is a grouping of trees and their fileglobs, and,...
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > fLoopManager
< The RLoopManager at the root of this computation graph. Never null.
RDataSource * fDataSource
Non-owning pointer to a data-source object. Null if no data-source. RLoopManager has ownership of the...
RDFDetail::RLoopManager * GetLoopManager() const
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
RDataFrame(std::string_view treeName, std::string_view filenameglob, const ColumnNames_t &defaultColumns={})
Build the dataframe.
ROOT::RDF::ColumnNames_t ColumnNames_t
Describe directory structure in memory.
A TTree represents a columnar dataset.
ROOT::RDataFrame FromSpec(const std::string &jsonFile)
Factory method to create an RDataFrame from a JSON specification file.
std::vector< std::string > ColumnNames_t
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