11#ifndef ROOT_RDF_TINTERFACE
12#define ROOT_RDF_TINTERFACE
37#include <initializer_list>
71template <
typename Proxied,
typename DataSource>
74using RNode = RInterface<::ROOT::Detail::RDF::RNodeBase, void>;
88template <
typename Proxied,
typename DataSource =
void>
98 template <
typename T,
typename W>
125 template <typename T = Proxied, typename std::enable_if<std::is_same<T, RLoopManager>::value,
int>
::type = 0>
185 template <typename F, typename std::enable_if<!std::is_convertible<F, std::string>::value,
int>
::type = 0>
189 RDFInternal::CheckFilter(
f);
190 using ColTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
191 constexpr auto nColumns = ColTypes_t::list_size;
193 const auto newColumns =
198 auto filterPtr = std::make_shared<F_t>(std::move(
f), validColumnNames,
fProxiedPtr, newColumns,
name);
211 template <typename F, typename std::enable_if<!std::is_convertible<F, std::string>::value,
int>
::type = 0>
227 template <
typename F>
252 using BaseNodeType_t =
typename std::remove_pointer<
decltype(upcastNodeOnHeap)>::type::element_type;
254 const auto jittedFilter = std::make_shared<RDFDetail::RJittedFilter>(
fLoopManager,
name);
294 template <typename F, typename std::enable_if<!std::is_convertible<F, std::string>::value,
int>
::type = 0>
297 return DefineImpl<F, RDFDetail::CustomColExtraArgs::None>(
name, std::move(expression), columns);
323 template <
typename F>
326 return DefineImpl<F, RDFDetail::CustomColExtraArgs::Slot>(
name, std::move(expression), columns);
353 template <
typename F>
356 return DefineImpl<F, RDFDetail::CustomColExtraArgs::SlotAndEntry>(
name, std::move(expression), columns);
378 auto jittedCustomColumn =
461 template <
typename... ColumnTypes>
466 return SnapshotImpl<ColumnTypes...>(treename, filename, columnList, options);
487 if (columnList.empty()) {
488 auto nEntries = *this->
Count();
489 auto snapshotRDF = std::make_shared<RInterface<RLoopManager>>(std::make_shared<RLoopManager>(nEntries));
494 std::stringstream snapCall;
502 snapCall <<
"*reinterpret_cast<ROOT::RDF::RResultPtr<ROOT::RDF::RInterface<ROOT::Detail::RDF::RLoopManager>>*>("
504 <<
") = reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RNodeBase>*>("
508 const auto dontConvertVector =
false;
512 for (
auto &
c : validColumnNames) {
513 const auto isCustom = std::find(customCols.begin(), customCols.end(),
c) != customCols.end();
519 if (!columnList.empty())
520 snapCall.seekp(-2, snapCall.cur);
521 snapCall <<
">(\"" << treename <<
"\", \"" << filename <<
"\", "
522 <<
"*reinterpret_cast<std::vector<std::string>*>("
553 return Snapshot(treename, filename, selectedColumns, options);
571 std::initializer_list<std::string> columnList,
575 return Snapshot(treename, filename, selectedColumns, options);
608 template <
typename... ColumnTypes>
611 auto staticSeq = std::make_index_sequence<
sizeof...(ColumnTypes)>();
612 return CacheImpl<ColumnTypes...>(columnList, staticSeq);
625 if (columnList.empty()) {
626 auto nEntries = *this->
Count();
633 std::stringstream cacheCall;
640 cacheCall <<
"*reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RLoopManager>*>("
642 <<
") = reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RNodeBase>*>("
646 for (
auto &
c : columnList) {
647 const auto isCustom = std::find(customCols.begin(), customCols.end(),
c) != customCols.end();
653 if (!columnList.empty())
654 cacheCall.seekp(-2, cacheCall.cur);
655 cacheCall <<
">(*reinterpret_cast<std::vector<std::string>*>("
674 columnNameRegexp,
"Cache");
675 return Cache(selectedColumns);
687 return Cache(selectedColumns);
711 if (stride == 0 || (end != 0 && end < begin))
712 throw std::runtime_error(
"Range: stride must be strictly greater than 0 and end must be greater than begin.");
716 auto rangePtr = std::make_shared<Range_t>(begin, end, stride,
fProxiedPtr);
749 template <
typename F>
752 using arg_types =
typename TTraits::CallableTraits<
decltype(
f)>::arg_types_nodecay;
753 using ret_type =
typename TTraits::CallableTraits<
decltype(
f)>::ret_type;
754 ForeachSlot(RDFInternal::AddSlotParameter<ret_type>(
f, arg_types()), columns);
779 template <
typename F>
782 using ColTypes_t = TypeTraits::RemoveFirstParameter_t<typename TTraits::CallableTraits<F>::arg_types>;
783 constexpr auto nColumns = ColTypes_t::list_size;
787 auto newColumns =
CheckAndFillDSColumns(validColumnNames, std::make_index_sequence<nColumns>(), ColTypes_t());
789 using Helper_t = RDFInternal::ForeachSlotHelper<F>;
793 std::make_unique<Action_t>(Helper_t(std::move(
f)), validColumnNames,
fProxiedPtr, std::move(newColumns));
829 template <typename F, typename T = typename TTraits::CallableTraits<F>::ret_type>
833 std::is_default_constructible<T>::value,
834 "reduce object cannot be default-constructed. Please provide an initialisation value (redIdentity)");
835 return Reduce(std::move(
f), columnName,
T());
852 template <typename F, typename T = typename TTraits::CallableTraits<F>::ret_type>
874 auto cSPtr = std::make_shared<ULong64_t>(0);
875 using Helper_t = RDFInternal::CountHelper;
903 template <
typename T,
typename COLL = std::vector<T>>
912 using Helper_t = RDFInternal::TakeHelper<T, T, COLL>;
914 auto valuesPtr = std::make_shared<COLL>();
918 std::make_unique<Action_t>(Helper_t(valuesPtr, nSlots), validColumnNames,
fProxiedPtr, std::move(newColumns));
945 template <
typename V = RDFDetail::RInferredType>
952 std::shared_ptr<::TH1D>
h(
nullptr);
954 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
955 h = model.GetHistogram();
956 h->SetDirectory(
nullptr);
959 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
960 RDFInternal::HistoUtils<::TH1D>::SetCanExtendAllAxes(*
h);
961 return CreateAction<RDFInternal::ActionTags::Histo1D, V>(validatedColumns,
h);
982 template <
typename V = RDFDetail::RInferredType>
985 const auto h_name = std::string(vName);
986 const auto h_title = h_name +
";" + h_name +
";count";
987 return Histo1D<V>({h_name.c_str(), h_title.c_str(), 128u, 0., 0.}, vName);
1009 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1012 const std::vector<std::string_view> columnViews = {vName, wName};
1016 std::shared_ptr<::TH1D>
h(
nullptr);
1018 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1021 return CreateAction<RDFInternal::ActionTags::Histo1D, V, W>(userColumns,
h);
1044 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1048 std::string str_vName{vName};
1049 std::string str_wName{wName};
1050 const auto h_name = str_vName +
"_weighted_" + str_wName;
1051 const auto h_title = str_vName +
", weights: " + str_wName +
";" + str_vName +
";count * " + str_wName;
1052 return Histo1D<V, W>({h_name.c_str(), h_title.c_str(), 128u, 0., 0.}, vName, wName);
1064 template <
typename V,
typename W>
1067 return Histo1D<V, W>(model,
"",
"");
1094 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
1097 std::shared_ptr<::TH2D>
h(
nullptr);
1099 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1102 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
1103 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
1105 const std::vector<std::string_view> columnViews = {v1Name, v2Name};
1109 return CreateAction<RDFInternal::ActionTags::Histo2D, V1, V2>(userColumns,
h);
1135 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1136 typename W = RDFDetail::RInferredType>
1140 std::shared_ptr<::TH2D>
h(
nullptr);
1142 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1145 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
1146 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
1148 const std::vector<std::string_view> columnViews = {v1Name, v2Name, wName};
1152 return CreateAction<RDFInternal::ActionTags::Histo2D, V1, V2, W>(userColumns,
h);
1155 template <
typename V1,
typename V2,
typename W>
1158 return Histo2D<V1, V2, W>(model,
"",
"",
"");
1185 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1186 typename V3 = RDFDetail::RInferredType>
1190 std::shared_ptr<::TH3D>
h(
nullptr);
1192 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1195 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
1196 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
1198 const std::vector<std::string_view> columnViews = {v1Name, v2Name, v3Name};
1202 return CreateAction<RDFInternal::ActionTags::Histo3D, V1, V2, V3>(userColumns,
h);
1232 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1233 typename V3 = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1237 std::shared_ptr<::TH3D>
h(
nullptr);
1239 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1242 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
1243 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
1245 const std::vector<std::string_view> columnViews = {v1Name, v2Name, v3Name, wName};
1249 return CreateAction<RDFInternal::ActionTags::Histo3D, V1, V2, V3, W>(userColumns,
h);
1252 template <
typename V1,
typename V2,
typename V3,
typename W>
1255 return Histo3D<V1, V2, V3, W>(model,
"",
"",
"",
"");
1283 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
1286 auto graph = std::make_shared<::TGraph>();
1287 const std::vector<std::string_view> columnViews = {v1Name, v2Name};
1295 if (!(validatedColumns[0].empty() && validatedColumns[1].empty())) {
1296 const auto g_name = std::string(v1Name) +
"_vs_" + std::string(v2Name);
1297 const auto g_title = std::string(v1Name) +
" vs " + std::string(v2Name);
1298 graph->SetNameTitle(g_name.c_str(), g_title.c_str());
1299 graph->GetXaxis()->SetTitle(std::string(v1Name).c_str());
1300 graph->GetYaxis()->SetTitle(std::string(v2Name).c_str());
1303 return CreateAction<RDFInternal::ActionTags::Graph, V1, V2>(validatedColumns,
graph);
1326 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
1330 std::shared_ptr<::TProfile>
h(
nullptr);
1332 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1336 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
1337 throw std::runtime_error(
"Profiles with no axes limits are not supported yet.");
1339 const std::vector<std::string_view> columnViews = {v1Name, v2Name};
1343 return CreateAction<RDFInternal::ActionTags::Profile1D, V1, V2>(userColumns,
h);
1369 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1370 typename W = RDFDetail::RInferredType>
1374 std::shared_ptr<::TProfile>
h(
nullptr);
1376 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1380 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
1381 throw std::runtime_error(
"Profile histograms with no axes limits are not supported yet.");
1383 const std::vector<std::string_view> columnViews = {v1Name, v2Name, wName};
1387 return CreateAction<RDFInternal::ActionTags::Profile1D, V1, V2, W>(userColumns,
h);
1390 template <
typename V1,
typename V2,
typename W>
1393 return Profile1D<V1, V2, W>(model,
"",
"",
"");
1420 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1421 typename V3 = RDFDetail::RInferredType>
1425 std::shared_ptr<::TProfile2D>
h(
nullptr);
1427 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1431 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
1432 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
1434 const std::vector<std::string_view> columnViews = {v1Name, v2Name, v3Name};
1438 return CreateAction<RDFInternal::ActionTags::Profile2D, V1, V2, V3>(userColumns,
h);
1467 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType,
1468 typename V3 = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1472 std::shared_ptr<::TProfile2D>
h(
nullptr);
1474 ROOT::Internal::RDF::RIgnoreErrorLevelRAII iel(
kError);
1478 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
1479 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
1481 const std::vector<std::string_view> columnViews = {v1Name, v2Name, v3Name, wName};
1485 return CreateAction<RDFInternal::ActionTags::Profile2D, V1, V2, V3, W>(userColumns,
h);
1488 template <
typename V1,
typename V2,
typename V3,
typename W>
1491 return Profile2D<V1, V2, V3, W>(model,
"",
"",
"",
"");
1518 template <
typename FirstColumn,
typename... OtherColumns,
typename T>
1521 auto h = std::make_shared<T>(std::forward<T>(model));
1522 if (!RDFInternal::HistoUtils<T>::HasAxisLimits(*
h)) {
1523 throw std::runtime_error(
"The absence of axes limits is not supported yet.");
1525 return CreateAction<RDFInternal::ActionTags::Fill, FirstColumn, OtherColumns...>(columnList,
h);
1547 template <
typename T>
1550 auto h = std::make_shared<T>(std::forward<T>(model));
1551 if (!RDFInternal::HistoUtils<T>::HasAxisLimits(*
h)) {
1552 throw std::runtime_error(
"The absence of axes limits is not supported yet.");
1554 return CreateAction<RDFInternal::ActionTags::Fill, RDFDetail::RInferredType>(bl,
h, bl.size());
1572 template<
typename V = RDFDetail::RInferredType>
1576 if (!value.empty()) {
1577 columns.emplace_back(std::string(value));
1580 if (std::is_same<V, RDFDetail::RInferredType>::value) {
1584 return Fill<V>(
TStatistic(), validColumnNames);
1605 template<
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1608 ColumnNames_t columns {std::string(value), std::string(weight)};
1609 constexpr auto vIsInferred = std::is_same<V, RDFDetail::RInferredType>::value;
1610 constexpr auto wIsInferred = std::is_same<W, RDFDetail::RInferredType>::value;
1616 if (vIsInferred && wIsInferred) {
1618 }
else if (vIsInferred != wIsInferred) {
1619 std::string error(
"The ");
1620 error += vIsInferred ?
"value " :
"weight ";
1621 error +=
"column type is explicit, while the ";
1622 error += vIsInferred ?
"weight " :
"value ";
1623 error +=
" is specified to be inferred. This case is not supported: please specify both types or none.";
1624 throw std::runtime_error(error);
1626 return Fill<V, W>(
TStatistic(), validColumnNames);
1651 template <
typename T = RDFDetail::RInferredType>
1655 using RetType_t = RDFDetail::MinReturnType_t<T>;
1656 auto minV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::max());
1657 return CreateAction<RDFInternal::ActionTags::Min, T>(userColumns, minV);
1681 template <
typename T = RDFDetail::RInferredType>
1685 using RetType_t = RDFDetail::MaxReturnType_t<T>;
1686 auto maxV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::lowest());
1687 return CreateAction<RDFInternal::ActionTags::Max, T>(userColumns, maxV);
1710 template <
typename T = RDFDetail::RInferredType>
1714 auto meanV = std::make_shared<double>(0);
1715 return CreateAction<RDFInternal::ActionTags::Mean, T>(userColumns, meanV);
1738 template <
typename T = RDFDetail::RInferredType>
1742 auto stdDeviationV = std::make_shared<double>(0);
1743 return CreateAction<RDFInternal::ActionTags::StdDev, T>(userColumns, stdDeviationV);
1769 template <
typename T = RDFDetail::RInferredType>
1772 const RDFDetail::SumReturnType_t<T> &initValue = RDFDetail::SumReturnType_t<T>{})
1775 auto sumV = std::make_shared<RDFDetail::SumReturnType_t<T>>(initValue);
1776 return CreateAction<RDFInternal::ActionTags::Sum, T>(userColumns, sumV);
1806 bool returnEmptyReport =
false;
1813 returnEmptyReport =
true;
1815 auto rep = std::make_shared<RCutFlowReport>();
1816 using Helper_t = RDFInternal::ReportHelper<Proxied>;
1845 allColumns.emplace_back(colName);
1850 std::for_each(columnNames.begin(), columnNames.end(), addIfNotInternal);
1855 allColumns.insert(allColumns.end(), branchNames.begin(), branchNames.end());
1860 allColumns.insert(allColumns.end(), dsColNames.begin(), dsColNames.end());
1882 const bool convertVector2RVec =
true;
1883 const auto isCustom = std::find(customCols.begin(), customCols.end(), column) != customCols.end();
1887 convertVector2RVec);
1891 const auto call =
"ROOT::Internal::RDF::TypeID2TypeName(typeid(__rdf" + std::to_string(
fLoopManager->
GetID()) +
1892 "::" + std::string(column) + colID +
"_type))";
1895 return *
reinterpret_cast<std::string *
>(calcRes);
1936 for (
auto column : columns) {
1938 definedColumns.emplace_back(column.first);
1941 return definedColumns;
1965 const auto branchNamesEnd = branchNames.end();
1966 if (branchNamesEnd != std::find(branchNames.begin(), branchNamesEnd, columnName))
2035 template <typename AccFun, typename MergeFun, typename R = typename TTraits::CallableTraits<AccFun>::ret_type,
2036 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
2037 typename ArgTypesNoDecay =
typename TTraits::CallableTraits<AccFun>::arg_types_nodecay,
2038 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
2039 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
2042 RDFInternal::CheckAggregate<R, MergeFun>(ArgTypesNoDecay());
2044 constexpr auto nColumns = ArgTypes::list_size;
2048 auto newColumns =
CheckAndFillDSColumns(validColumnNames, std::make_index_sequence<nColumns>(), ArgTypes());
2050 auto accObjPtr = std::make_shared<U>(aggIdentity);
2051 using Helper_t = RDFInternal::AggregateHelper<AccFun, MergeFun, R, T, U>;
2053 auto action = std::make_unique<Action_t>(
2054 Helper_t(std::move(aggregator), std::move(merger), accObjPtr,
fLoopManager->
GetNSlots()), validColumnNames,
2073 template <typename AccFun, typename MergeFun, typename R = typename TTraits::CallableTraits<AccFun>::ret_type,
2074 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
2075 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
2076 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
2080 std::is_default_constructible<U>::value,
2081 "aggregated object cannot be default-constructed. Please provide an initialisation value (aggIdentity)");
2082 return Aggregate(std::move(aggregator), std::move(merger), columnName, U());
2119 template <
typename... ColumnTypes,
typename Helper>
2122 constexpr auto nColumns =
sizeof...(ColumnTypes);
2128 using AH = RDFDetail::RActionImpl<Helper>;
2129 static_assert(std::is_base_of<AH, Helper>::value && std::is_convertible<Helper *, AH *>::value,
2130 "Action helper of type T must publicly inherit from ROOT::Detail::RDF::RActionImpl<T>");
2133 auto resPtr = helper.GetResultPtr();
2138 auto action = std::make_unique<Action_t>(Helper(std::forward<Helper>(helper)), validColumnNames,
fProxiedPtr,
2167 template <
typename... ColumnTypes>
2172 auto displayer = std::make_shared<RDFInternal::RDisplay>(columnList,
GetColumnTypeNamesList(columnList), nRows);
2173 return CreateAction<RDFInternal::ActionTags::Display, ColumnTypes...>(columnList, displayer);
2187 auto displayer = std::make_shared<RDFInternal::RDisplay>(columnList,
GetColumnTypeNamesList(columnList), nRows);
2188 return CreateAction<RDFInternal::ActionTags::Display, RDFDetail::RInferredType>(columnList, displayer,
2204 columnNameRegexp,
"Display");
2205 return Display(selectedColumns, nRows);
2218 return Display(selectedColumns, nRows);
2227 const auto entryColName =
"rdfentry_";
2228 auto entryColGen = [](
unsigned int,
ULong64_t entry) {
return entry; };
2229 using NewColEntry_t =
2232 auto entryColumn = std::make_shared<NewColEntry_t>(
fLoopManager, entryColName, std::move(entryColGen),
2234 newCols.
AddName(entryColName);
2235 newCols.
AddColumn(entryColumn, entryColName);
2238 auto retTypeDeclaration =
"namespace __rdf" + std::to_string(
fLoopManager->
GetID()) +
" { using " + entryColName +
2239 std::to_string(entryColumn->
GetID()) +
"_type = ULong64_t; }";
2243 const auto slotColName =
"rdfslot_";
2244 auto slotColGen = [](
unsigned int slot) {
return slot; };
2247 auto slotColumn = std::make_shared<NewColSlot_t>(
fLoopManager, slotColName, std::move(slotColGen),
2251 newCols.
AddColumn(slotColumn, slotColName);
2256 retTypeDeclaration =
"namespace __rdf" + std::to_string(
fLoopManager->
GetID()) +
" { using " + slotColName +
2257 std::to_string(slotColumn->
GetID()) +
"_type = unsigned int; }";
2268 std::vector<std::string> types;
2270 for (
auto column : columnList) {
2279 std::string error(callerName);
2280 error +=
" was called with ImplicitMT enabled, but multi-thread is not supported.";
2281 throw std::runtime_error(error);
2286 template <
typename ActionTag,
typename... BranchTypes,
typename ActionResultType,
2287 typename std::enable_if<!RDFInternal::TNeedJitting<BranchTypes...>::value,
int>
::type = 0>
2290 constexpr auto nColumns =
sizeof...(BranchTypes);
2299 auto action = RDFInternal::BuildAction<BranchTypes...>(validColumnNames,
r, nSlots,
fProxiedPtr, ActionTag{},
2300 std::move(newColumns));
2308 template <
typename ActionTag,
typename... BranchTypes,
typename ActionResultType,
2309 typename std::enable_if<RDFInternal::TNeedJitting<BranchTypes...>::value,
int>
::type = 0>
2313 auto realNColumns = (nColumns > -1 ? nColumns :
sizeof...(BranchTypes));
2319 auto rOnHeap = RDFInternal::MakeSharedOnHeap(
r);
2322 using BaseNodeType_t =
typename std::remove_pointer<
decltype(upcastNodeOnHeap)>::type::element_type;
2325 auto jittedActionOnHeap =
2326 RDFInternal::MakeSharedOnHeap(std::make_shared<RDFInternal::RJittedAction>(*
fLoopManager));
2329 validColumnNames, upcastNodeOnHeap,
typeid(std::shared_ptr<ActionResultType>),
typeid(ActionTag), rOnHeap,
2336 template <typename F, typename CustomColumnType, typename RetType = typename TTraits::CallableTraits<F>::ret_type>
2344 using ArgTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
2345 using ColTypesTmp_t =
typename RDFInternal::RemoveFirstParameterIf<
2346 std::is_same<CustomColumnType, RDFDetail::CustomColExtraArgs::Slot>::value, ArgTypes_t>
::type;
2347 using ColTypes_t =
typename RDFInternal::RemoveFirstTwoParametersIf<
2348 std::is_same<CustomColumnType, RDFDetail::CustomColExtraArgs::SlotAndEntry>::value, ColTypesTmp_t>
::type;
2350 constexpr auto nColumns = ColTypes_t::list_size;
2354 auto newColumns =
CheckAndFillDSColumns(validColumnNames, std::make_index_sequence<nColumns>(), ColTypes_t());
2358 auto newColumn = std::make_shared<NewCol_t>(
fLoopManager,
name, std::forward<F>(expression), validColumnNames,
2363 if (retTypeName.empty()) {
2368 retTypeName =
"void /* The type of column \"" + std::string(
name) +
"\" (" + demangledType +
2369 ") is not known to the interpreter. */";
2371 const auto retTypeDeclaration =
"namespace __rdf" + std::to_string(
fLoopManager->
GetID()) +
2372 " { " + +
" using " + std::string(
name) + std::to_string(newColumn->GetID()) +
2373 "_type = " + retTypeName +
"; }";
2381 return newInterface;
2387 template <typename F, typename CustomColumnType, typename RetType = typename TTraits::CallableTraits<F>::ret_type,
2388 bool IsFStringConv = std::is_convertible<F, std::string>::value,
2389 bool IsRetTypeDefConstr = std::is_default_constructible<RetType>::value>
2390 typename std::enable_if<!IsFStringConv && !IsRetTypeDefConstr, RInterface<Proxied, DS_t>>
::type
2393 static_assert(std::is_default_constructible<typename TTraits::CallableTraits<F>::ret_type>::value,
2394 "Error in `Define`: type returned by expression is not default-constructible");
2408 template <
typename... ColumnTypes>
2419 const std::string fullTreename(treename);
2421 const auto lastSlash = treename.rfind(
'/');
2423 if (std::string_view::npos != lastSlash) {
2424 dirname = treename.substr(0, lastSlash);
2425 treename = treename.substr(lastSlash + 1, treename.size());
2429 std::unique_ptr<RDFInternal::RActionBase> actionPtr;
2434 actionPtr.reset(
new Action_t(Helper_t(filename, dirname, treename, validCols, columnList, options), validCols,
2440 actionPtr.reset(
new Action_t(
2441 Helper_t(
fLoopManager->
GetNSlots(), filename, dirname, treename, validCols, columnList, options), validCols,
2448 std::move(actionPtr));
2453 template <
typename... BranchTypes, std::size_t...
S>
2457 constexpr bool areCopyConstructible =
2458 RDFInternal::TEvalAnd<std::is_copy_constructible<BranchTypes>::value...>::value;
2459 static_assert(areCopyConstructible,
"Columns of a type which is not copy constructible cannot be cached yet.");
2465 auto colHolders = std::make_tuple(Take<BranchTypes>(columnList[
S])...);
2466 auto ds = std::make_unique<
RLazyDS<BranchTypes...>>(std::make_pair(columnList[
S], std::get<S>(colHolders))...);
2494 template <
typename... ColumnTypes, std::size_t...
S>
unsigned long long ULong64_t
typedef void((*Func_t)())
unsigned int GetID() const
Return the unique identifier of this RCustomColumnBase.
The head node of a RDF computation graph.
const std::map< std::string, std::string > & GetAliasMap() const
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
void ToJitDeclare(const std::string &s)
void ToJitExec(const std::string &s)
void Run()
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
unsigned int GetID() const
void JitDeclarations()
Declare to the interpreter type aliases and other entities required by RDF jitted nodes.
void AddColumnAlias(const std::string &alias, const std::string &colName)
RDataSource * GetDataSource() const
unsigned int GetNSlots() const
void Book(RDFInternal::RActionBase *actionPtr)
Helper class that provides the operation graph nodes.
An action node in a RDF computation graph.
Encapsulates the columns defined by the user.
void AddColumn(const std::shared_ptr< RDFDetail::RCustomColumnBase > &column, std::string_view name)
Internally it recreates the map with the new column, and swaps with the old one.
ColumnNames_t GetNames() const
Returns the list of the names of the defined columns.
bool HasName(std::string_view name) const
Check if the provided name is tracked in the names list.
const RCustomColumnBasePtrMap_t & GetColumns() const
Returns the list of the pointers to the defined columns.
void AddName(std::string_view name)
Internally it recreates the map with the new column name, and swaps with the old one.
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
virtual const std::vector< std::string > & GetColumnNames() const =0
Returns a reference to the collection of the dataset's column names.
virtual bool HasColumn(std::string_view) const =0
Checks if the dataset has a certain column.
The public interface to the RDataFrame federation of classes.
RInterface(const RInterface &)=default
Copy-ctor for RInterface.
RResultPtr<::TH1D > Histo1D(std::string_view vName, std::string_view wName)
Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
RResultPtr<::TH1D > Histo1D(const TH1DModel &model={"", "", 128u, 0., 0.})
Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
RResultPtr<::TH2D > Histo2D(const TH2DModel &model)
RResultPtr<::TProfile > Profile1D(const TProfile1DModel &model, std::string_view v1Name="", std::string_view v2Name="")
Fill and return a one-dimensional profile (lazy action)
RInterface(const std::shared_ptr< Proxied > &proxied, RLoopManager &lm, const RDFInternal::RBookedCustomColumns &columns, RDataSource *ds)
RResultPtr< RInterface< RLoopManager > > Snapshot(std::string_view treename, std::string_view filename, std::string_view columnNameRegexp="", const RSnapshotOptions &options=RSnapshotOptions())
Save selected columns to disk, in a new TTree treename in file filename.
RResultPtr< TStatistic > Stats(std::string_view value="")
Return a TStatistic object, filled once per event (lazy action)
RLoopManager * GetLoopManager() const
RResultPtr<::TGraph > Graph(std::string_view v1Name="", std::string_view v2Name="")
Fill and return a graph (lazy action)
RInterface< Proxied, DS_t > DefineSlot(std::string_view name, F expression, const ColumnNames_t &columns={})
Creates a custom column with a value dependent on the processing slot.
RResultPtr< double > StdDev(std::string_view columnName="")
Return the unbiased standard deviation of processed column values (lazy action)
unsigned int GetNSlots() const
Gets the number of data processing slots.
RInterface(const std::shared_ptr< Proxied > &proxied)
Only enabled when building a RInterface<RLoopManager>
RResultPtr< T > Fill(T &&model, const ColumnNames_t &bl)
Return an object of type T on which T::Fill will be called once per event (lazy action)
RResultPtr< ActionResultType > CreateAction(const ColumnNames_t &columns, const std::shared_ptr< ActionResultType > &r, const int nColumns=-1)
void ForeachSlot(F f, const ColumnNames_t &columns={})
Execute a user-defined function requiring a processing slot index on each entry (instant action)
RInterface< RLoopManager > Cache(const ColumnNames_t &columnList)
Save selected columns in memory.
RResultPtr< RDisplay > Display(std::initializer_list< std::string > columnList, const int &nRows=5)
Provides a representation of the columns in the dataset.
RInterface< Proxied, DS_t > Define(std::string_view name, F expression, const ColumnNames_t &columns={})
Creates a custom column.
RResultPtr< TStatistic > Stats(std::string_view value, std::string_view weight)
Return a TStatistic object, filled once per event (lazy action)
RDataSource * fDataSource
Non-owning pointer to a data-source object. Null if no data-source. RLoopManager has ownership of the...
RResultPtr<::TH2D > Histo2D(const TH2DModel &model, std::string_view v1Name="", std::string_view v2Name="")
Fill and return a two-dimensional histogram (lazy action)
ColumnNames_t GetValidatedColumnNames(const unsigned int nColumns, const ColumnNames_t &columns)
Prepare the call to the GetValidatedColumnNames routine, making sure that GetBranchNames,...
RResultPtr<::TProfile > Profile1D(const TProfile1DModel &model)
RDFInternal::RBookedCustomColumns fCustomColumns
Contains the custom columns defined up to this node.
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, const std::initializer_list< std::string > &columns)
Append a filter to the call graph.
RResultPtr< double > Mean(std::string_view columnName="")
Return the mean of processed column values (lazy action)
RResultPtr< RInterface< RLoopManager > > Snapshot(std::string_view treename, std::string_view filename, std::initializer_list< std::string > columnList, const RSnapshotOptions &options=RSnapshotOptions())
Save selected columns to disk, in a new TTree treename in file filename.
std::enable_if<!IsFStringConv &&!IsRetTypeDefConstr, RInterface< Proxied, DS_t > >::type DefineImpl(std::string_view, F, const ColumnNames_t &)
RInterface< Proxied, DS_t > Alias(std::string_view alias, std::string_view columnName)
Allow to refer to a column with a different name.
RInterface< RLoopManager > Cache(const ColumnNames_t &columnList)
Save selected columns in memory.
RInterface< RLoopManager > Cache(std::string_view columnNameRegexp="")
Save selected columns in memory.
RResultPtr< RDisplay > Display(std::string_view columnNameRegexp="", const int &nRows=5)
Provides a representation of the columns in the dataset.
RLoopManager * fLoopManager
friend class RDFInternal::GraphDrawing::GraphCreatorHelper
RResultPtr<::TH2D > Histo2D(const TH2DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view wName)
Fill and return a weighted two-dimensional histogram (lazy action)
RInterface & operator=(const RInterface &)=default
Copy-assignment operator for RInterface.
RResultPtr< RDFDetail::SumReturnType_t< T > > Sum(std::string_view columnName="", const RDFDetail::SumReturnType_t< T > &initValue=RDFDetail::SumReturnType_t< T >{})
Return the sum of processed column values (lazy action)
RResultPtr< ULong64_t > Count()
Return the number of entries processed (lazy action)
RResultPtr< T > Fill(T &&model, const ColumnNames_t &columnList)
Return an object of type T on which T::Fill will be called once per event (lazy action)
RInterface< Proxied, DS_t > Define(std::string_view name, std::string_view expression)
Creates a custom column.
std::shared_ptr< Proxied > fProxiedPtr
Smart pointer to the graph node encapsulated by this RInterface.
RResultPtr<::TH1D > Histo1D(std::string_view vName)
Fill and return a one-dimensional histogram with the values of a column (lazy action)
RDFInternal::RBookedCustomColumns CheckAndFillDSColumns(ColumnNames_t validCols, std::index_sequence< S... >, TTraits::TypeList< ColumnTypes... >)
ColumnNames_t GetColumnNames()
Returns the names of the available columns.
RResultPtr< RInterface< RLoopManager > > SnapshotImpl(std::string_view treename, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options)
Implementation of snapshot.
RResultPtr< RDisplay > Display(const ColumnNames_t &columnList, const int &nRows=5)
Provides a representation of the columns in the dataset.
RResultPtr<::TH1D > Histo1D(const TH1DModel &model, std::string_view vName, std::string_view wName)
Fill and return a one-dimensional histogram with the weighted values of a column (lazy action)
RInterface< RDFDetail::RRange< Proxied >, DS_t > Range(unsigned int end)
Creates a node that filters entries based on range.
RResultPtr< COLL > Take(std::string_view column="")
Return a collection of values of a column (lazy action, returns a std::vector by default)
RInterface< RLoopManager > Cache(std::initializer_list< std::string > columnList)
Save selected columns in memory.
RResultPtr<::TProfile2D > Profile2D(const TProfile2DModel &model, std::string_view v1Name="", std::string_view v2Name="", std::string_view v3Name="")
Fill and return a two-dimensional profile (lazy action)
RResultPtr< typename Helper::Result_t > Book(Helper &&helper, const ColumnNames_t &columns={})
Book execution of a custom action using a user-defined helper object.
RResultPtr< RDisplay > Display(const ColumnNames_t &columnList, const int &nRows=5)
Provides a representation of the columns in the dataset.
const std::shared_ptr< Proxied > & GetProxiedPtr() const
RResultPtr<::TH3D > Histo3D(const TH3DModel &model, std::string_view v1Name="", std::string_view v2Name="", std::string_view v3Name="")
Fill and return a three-dimensional histogram (lazy action)
RResultPtr< ActionResultType > CreateAction(const ColumnNames_t &columns, const std::shared_ptr< ActionResultType > &r)
RResultPtr< RInterface< RLoopManager > > Snapshot(std::string_view treename, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options=RSnapshotOptions())
Save selected columns to disk, in a new TTree treename in file filename.
RResultPtr< RCutFlowReport > Report()
Gather filtering statistics.
RResultPtr<::TProfile2D > Profile2D(const TProfile2DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view v3Name, std::string_view wName)
Fill and return a two-dimensional profile (lazy action)
RResultPtr< RInterface< RLoopManager > > Snapshot(std::string_view treename, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options=RSnapshotOptions())
Save selected columns to disk, in a new TTree treename in file filename.
RResultPtr< U > Aggregate(AccFun aggregator, MergeFun merger, std::string_view columnName="")
Execute a user-defined accumulation operation on the processed column values in each processing slot.
RInterface< Proxied, DS_t > DefineSlotEntry(std::string_view name, F expression, const ColumnNames_t &columns={})
Creates a custom column with a value dependent on the processing slot and the current entry.
RResultPtr< RDFDetail::MinReturnType_t< T > > Min(std::string_view columnName="")
Return the minimum of processed column values (lazy action)
RResultPtr< T > Reduce(F f, std::string_view columnName="")
Execute a user-defined reduce operation on the values of a column.
void Foreach(F f, const ColumnNames_t &columns={})
Execute a user-defined function on each entry (instant action)
RInterface< RDFDetail::RJittedFilter, DS_t > Filter(std::string_view expression, std::string_view name="")
Append a filter to the call graph.
RResultPtr<::TProfile2D > Profile2D(const TProfile2DModel &model)
std::enable_if< std::is_default_constructible< RetType >::value, RInterface< Proxied, DS_t > >::type DefineImpl(std::string_view name, F &&expression, const ColumnNames_t &columns)
std::string GetColumnType(std::string_view column)
Return the type of a given column as a string.
ColumnNames_t GetDefinedColumnNames()
Returns the names of the defined columns.
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, const ColumnNames_t &columns={}, std::string_view name="")
Append a filter to the call graph.
RResultPtr< U > Aggregate(AccFun aggregator, MergeFun merger, std::string_view columnName, const U &aggIdentity)
Execute a user-defined accumulation operation on the processed column values in each processing slot.
RInterface(RInterface &&)=default
Move-ctor for RInterface.
RResultPtr< T > Reduce(F f, std::string_view columnName, const T &redIdentity)
Execute a user-defined reduce operation on the values of a column.
void CheckIMTDisabled(std::string_view callerName)
RResultPtr<::TH3D > Histo3D(const TH3DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view v3Name, std::string_view wName)
Fill and return a three-dimensional histogram (lazy action)
bool HasColumn(std::string_view columnName)
Checks if a column is present in the dataset.
RInterface< RLoopManager > CacheImpl(const ColumnNames_t &columnList, std::index_sequence< S... > s)
Implementation of cache.
RDFDetail::ColumnNames_t ColumnNames_t
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, std::string_view name)
Append a filter to the call graph.
RInterface< RDFDetail::RRange< Proxied >, DS_t > Range(unsigned int begin, unsigned int end, unsigned int stride=1)
Creates a node that filters entries based on range: [begin, end)
std::vector< std::string > GetColumnTypeNamesList(const ColumnNames_t &columnList)
std::vector< std::string > GetFilterNames()
Returns the names of the filters created.
RResultPtr<::TH1D > Histo1D(const TH1DModel &model={"", "", 128u, 0., 0.}, std::string_view vName="")
Fill and return a one-dimensional histogram with the values of a column (lazy action)
RResultPtr<::TProfile > Profile1D(const TProfile1DModel &model, std::string_view v1Name, std::string_view v2Name, std::string_view wName)
Fill and return a one-dimensional profile (lazy action)
RResultPtr<::TH3D > Histo3D(const TH3DModel &model)
RResultPtr< RDFDetail::MaxReturnType_t< T > > Max(std::string_view columnName="")
Return the maximum of processed column values (lazy action)
A RDataSource implementation which is built on top of result proxies.
Smart pointer for the return type of actions.
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
A Graph is a graphics object made of two arrays X and Y with npoints each.
Statistical variable, defined by its mean and variance (RMS).
basic_string_view< char > string_view
RResultPtr< T > MakeResultPtr(const std::shared_ptr< T > &r, RLoopManager &df, std::shared_ptr< ROOT::Internal::RDF::RActionBase > actionPtr)
Create a RResultPtr and set its pointer to the corresponding RAction This overload is invoked by non-...
ColumnNames_t GetBranchNames(TTree &t, bool allowDuplicates=true)
Get all the branches names, including the ones of the friend trees.
std::shared_ptr< RNodeBase > UpcastNode(std::shared_ptr< RNodeBase > ptr)
std::string TypeID2TypeName(const std::type_info &id)
Returns the name of a type starting from its type_info An empty string is returned in case of failure...
std::vector< std::string > GetFilterNames(const std::shared_ptr< RLoopManager > &loopManager)
ColumnNames_t ConvertRegexToColumns(const RDFInternal::RBookedCustomColumns &customColumns, TTree *tree, ROOT::RDF::RDataSource *dataSource, std::string_view columnNameRegexp, std::string_view callerName)
std::string PrettyPrintAddr(const void *const addr)
void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames)
bool AtLeastOneEmptyString(const std::vector< std::string_view > strings)
HeadNode_t CreateSnapshotRDF(const ColumnNames_t &validCols, std::string_view treeName, std::string_view fileName, bool isLazy, RLoopManager &loopManager, std::unique_ptr< RDFInternal::RActionBase > actionPtr)
std::string JitBuildAction(const ColumnNames_t &bl, void *prevNode, const std::type_info &art, const std::type_info &at, void *rOnHeap, TTree *tree, const unsigned int nSlots, const RDFInternal::RBookedCustomColumns &customCols, RDataSource *ds, std::shared_ptr< RJittedAction > *jittedActionOnHeap, unsigned int namespaceID)
bool IsInternalColumn(std::string_view colName)
Long64_t InterpreterCalc(const std::string &code, const std::string &context)
std::string ColumnName2ColumnTypeName(const std::string &colName, unsigned int namespaceID, TTree *tree, RDataSource *ds, bool isCustomColumn, bool vector2rvec, unsigned int customColID)
Return a string containing the type of the given branch.
ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns, const ColumnNames_t &validCustomColumns, RDataSource *ds)
Given the desired number of columns and the user-provided list of columns:
void BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm, RDataSource *ds, const std::shared_ptr< RJittedCustomColumn > &jittedCustomColumn, const RDFInternal::RBookedCustomColumns &customCols, const ColumnNames_t &branches)
void BookFilterJit(RJittedFilter *jittedFilter, void *prevNodeOnHeap, std::string_view name, std::string_view expression, const std::map< std::string, std::string > &aliasMap, const ColumnNames_t &branches, const RDFInternal::RBookedCustomColumns &customCols, TTree *tree, RDataSource *ds, unsigned int namespaceID)
void CheckCustomColumn(std::string_view definedCol, TTree *treePtr, const ColumnNames_t &customCols, const std::map< std::string, std::string > &aliasMap, const ColumnNames_t &dataSourceColumns)
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
ROOT type_traits extensions.
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
ROOT::Detail::RDF::ColumnNames_t ColumnNames_t
void DisableImplicitMT()
Disables the implicit multi-threading in ROOT (see EnableImplicitMT).
std::pair< Double_t, Double_t > Range_t
RooArgSet S(const RooAbsArg &v1)
char * DemangleTypeIdName(const std::type_info &ti, int &errorCode)
Demangle in a portable way the type id name.
static constexpr double s
A collection of options to steer the creation of the dataset on file.
bool fLazy
Delay the snapshot of the dataset.
A struct which stores the parameters of a TH1D.
std::shared_ptr<::TH1D > GetHistogram() const
A struct which stores the parameters of a TH2D.
std::shared_ptr<::TH2D > GetHistogram() const
A struct which stores the parameters of a TH3D.
std::shared_ptr<::TH3D > GetHistogram() const
A struct which stores the parameters of a TProfile.
std::shared_ptr<::TProfile > GetProfile() const
A struct which stores the parameters of a TProfile2D.
std::shared_ptr<::TProfile2D > GetProfile() const
Lightweight storage for a collection of types.