11#ifndef ROOT_RDF_TINTERFACE
12#define ROOT_RDF_TINTERFACE
49#include <initializer_list>
59#include <unordered_set>
83template <
typename Proxied,
typename DataSource>
115template <
typename Proxied,
typename DataSource =
void>
124 template <
typename T,
typename W>
226 RDFInternal::CheckFilter(
f);
227 using ColTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
228 constexpr auto nColumns = ColTypes_t::list_size;
262 template <
typename F>
348 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
399 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
474 template <
typename F>
504 template <
typename F>
535 constexpr auto where =
"Define";
583 template <
typename F>
602 template <
typename F>
606 "RedefineSlotEntry");
625 constexpr auto where =
"Redefine";
676 template <
typename T>
679 constexpr auto where{
"DefaultValueFor"};
697 auto newColumn = std::make_shared<ROOT::Internal::RDF::RDefaultValueFor<T>>(
866 template <
typename F>
905 template <
typename F>
953 template <
typename F>
975 template <
typename F>
1015 template <
typename F>
1046 template <
typename F>
1228 constexpr auto where =
"Alias";
1349 throw std::runtime_error(
"Snapshotting from TTree to RNTuple is not yet supported. The current recommended "
1350 "way to convert TTrees to RNTuple is through the RNTupleImporter.");
1358 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1369 throw std::runtime_error(
1370 "RDataFrame: Cannot snapshot to RNTuple - this installation of ROOT has not been build with ROOT7 "
1371 "components enabled.");
1377 "The default Snapshot output data format is TTree, but the input data format is RNTuple. If you "
1378 "want to Snapshot to RNTuple or suppress this warning, set the appropriate fOutputFormat option in "
1379 "RSnapshotOptions. Note that this current default behaviour might change in the future.");
1384 auto newRDF = std::make_shared<RInterface<RLoopManager>>(
1387 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1426 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1461 std::initializer_list<std::string>
columnList,
1506 auto staticSeq = std::make_index_sequence<
sizeof...(ColumnTypes)>();
1533 cacheCall <<
"*reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RLoopManager>*>("
1535 <<
") = reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RNodeBase>*>("
1548 cacheCall <<
">(*reinterpret_cast<std::vector<std::string>*>("
1575 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1618 if (
stride == 0 || (end != 0 && end < begin))
1619 throw std::runtime_error(
"Range: stride must be strictly greater than 0 and end must be greater than begin.");
1655 template <
typename F>
1658 using arg_types =
typename TTraits::CallableTraits<
decltype(
f)>::arg_types_nodecay;
1659 using ret_type =
typename TTraits::CallableTraits<
decltype(
f)>::ret_type;
1685 template <
typename F>
1689 constexpr auto nColumns = ColTypes_t::list_size;
1694 using Helper_t = RDFInternal::ForeachSlotHelper<F>;
1736 std::is_default_constructible<T>::value,
1737 "reduce object cannot be default-constructed. Please provide an initialisation value (redIdentity)");
1777 auto cSPtr = std::make_shared<ULong64_t>(0);
1778 using Helper_t = RDFInternal::CountHelper;
1805 template <
typename T,
typename COLL = std::vector<T>>
1813 using Helper_t = RDFInternal::TakeHelper<T, T, COLL>;
1815 auto valuesPtr = std::make_shared<COLL>();
1848 template <
typename V = RDFDetail::RInferredType>
1855 std::shared_ptr<::TH1D>
h(
nullptr);
1857 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
1858 h = model.GetHistogram();
1859 h->SetDirectory(
nullptr);
1862 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
1863 RDFInternal::HistoUtils<::TH1D>::SetCanExtendAllAxes(*
h);
1884 template <
typename V = RDFDetail::RInferredType>
1910 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1917 std::shared_ptr<::TH1D>
h(
nullptr);
1919 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
1944 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1964 template <
typename V,
typename W>
1998 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2001 std::shared_ptr<::TH2D>
h(
nullptr);
2003 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2006 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2007 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2044 std::shared_ptr<::TH2D>
h(
nullptr);
2046 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2049 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2050 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2059 template <
typename V1,
typename V2,
typename W>
2095 std::string_view
v3Name =
"")
2097 std::shared_ptr<::TH3D>
h(
nullptr);
2099 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2102 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2103 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2146 std::shared_ptr<::TH3D>
h(
nullptr);
2148 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2151 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2152 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2161 template <
typename V1,
typename V2,
typename V3,
typename W>
2192 std::shared_ptr<::THnD>
h(
nullptr);
2194 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2197 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2199 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2200 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2226 std::shared_ptr<::THnD>
h(
nullptr);
2228 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2231 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2233 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2234 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2269 template <
typename X = RDFDetail::RInferredType,
typename Y = RDFDetail::RInferredType>
2272 auto graph = std::make_shared<::TGraph>();
2336 std::string_view
exh =
"", std::string_view
eyl =
"", std::string_view
eyh =
"")
2338 auto graph = std::make_shared<::TGraphAsymmErrors>();
2380 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2384 std::shared_ptr<::TProfile>
h(
nullptr);
2386 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2390 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2391 throw std::runtime_error(
"Profiles with no axes limits are not supported yet.");
2429 std::shared_ptr<::TProfile>
h(
nullptr);
2431 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2435 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2436 throw std::runtime_error(
"Profile histograms with no axes limits are not supported yet.");
2448 template <
typename V1,
typename V2,
typename W>
2484 std::string_view
v2Name =
"", std::string_view
v3Name =
"")
2486 std::shared_ptr<::TProfile2D>
h(
nullptr);
2488 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2492 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2493 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2534 std::shared_ptr<::TProfile2D>
h(
nullptr);
2536 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2540 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2541 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2552 template <
typename V1,
typename V2,
typename V3,
typename W>
2595 auto h = std::make_shared<std::decay_t<T>>(std::forward<T>(model));
2596 if (!RDFInternal::HistoUtils<T>::HasAxisLimits(*
h)) {
2597 throw std::runtime_error(
"The absence of axes limits is not supported yet.");
2618 template <
typename V = RDFDetail::RInferredType>
2622 if (!
value.empty()) {
2626 if (std::is_same<V, RDFDetail::RInferredType>::value) {
2650 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2654 constexpr auto vIsInferred = std::is_same<V, RDFDetail::RInferredType>::value;
2655 constexpr auto wIsInferred = std::is_same<W, RDFDetail::RInferredType>::value;
2664 std::string error(
"The ");
2666 error +=
"column type is explicit, while the ";
2668 error +=
" is specified to be inferred. This case is not supported: please specify both types or none.";
2669 throw std::runtime_error(error);
2696 template <
typename T = RDFDetail::RInferredType>
2700 using RetType_t = RDFDetail::MinReturnType_t<T>;
2701 auto minV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::max());
2726 template <
typename T = RDFDetail::RInferredType>
2730 using RetType_t = RDFDetail::MaxReturnType_t<T>;
2731 auto maxV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::lowest());
2755 template <
typename T = RDFDetail::RInferredType>
2759 auto meanV = std::make_shared<double>(0);
2783 template <
typename T = RDFDetail::RInferredType>
2814 template <
typename T = RDFDetail::RInferredType>
2817 const RDFDetail::SumReturnType_t<T> &
initValue = RDFDetail::SumReturnType_t<T>{})
2820 auto sumV = std::make_shared<RDFDetail::SumReturnType_t<T>>(
initValue);
2860 auto rep = std::make_shared<RCutFlowReport>();
2861 using Helper_t = RDFInternal::ReportHelper<Proxied>;
2932 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
2933 typename ArgTypesNoDecay =
typename TTraits::CallableTraits<AccFun>::arg_types_nodecay,
2934 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
2935 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
2945 using Helper_t = RDFInternal::AggregateHelper<AccFun, MergeFun, R, T, U>;
2947 auto action = std::make_unique<Action_t>(
2967 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
2968 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
2969 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
2973 std::is_default_constructible<U>::value,
2974 "aggregated object cannot be default-constructed. Please provide an initialisation value (aggIdentity)");
3044 using HelperT = std::decay_t<Helper>;
3047 static_assert(std::is_base_of<AH, HelperT>::value && std::is_convertible<HelperT *, AH *>::value,
3048 "Action helper of type T must publicly inherit from ROOT::Detail::RDF::RActionImpl<T>");
3050 auto hPtr = std::make_shared<HelperT>(std::forward<Helper>(
helper));
3053 if (std::is_same<FirstColumn, RDFDetail::RInferredType>::value &&
columns.empty()) {
3160 if (
where.compare(0, 8,
"Redefine") != 0) {
3170 using ArgTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
3172 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::Slot>::value,
ArgTypes_t>
::type;
3174 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::SlotAndEntry>::value,
ColTypesTmp_t>
::type;
3176 constexpr auto nColumns = ColTypes_t::list_size;
3206 bool IsFStringConv = std::is_convertible<F, std::string>::value,
3208 std::enable_if_t<!IsFStringConv && !IsRetTypeDefConstr, RInterface<Proxied, DS_t>>
3211 static_assert(std::is_default_constructible<typename TTraits::CallableTraits<F>::ret_type>
::value,
3212 "Error in `Define`: type returned by expression is not default-constructible");
3239 throw std::runtime_error(
"Snapshotting from TTree to RNTuple is not yet supported. The current recommended "
3240 "way to convert TTrees to RNTuple is through the RNTupleImporter.");
3246 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
3256 throw std::runtime_error(
3257 "RDataFrame: Cannot snapshot to RNTuple - this installation of ROOT has not been build with ROOT7 "
3258 "components enabled.");
3264 "The default Snapshot output data format is TTree, but the input data format is RNTuple. If you "
3265 "want to Snapshot to RNTuple or suppress this warning, set the appropriate fOutputFormat option in "
3266 "RSnapshotOptions. Note that this current default behaviour might change in the future.");
3274 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
3292 template <
typename...
ColTypes, std::size_t... S>
3299 RDFInternal::TEvalAnd<std::is_copy_constructible<ColTypes>::value...>
::value;
3300 static_assert(
areCopyConstructible,
"Columns of a type which is not copy constructible cannot be cached yet.");
3313 template <
bool IsSingleColumn,
typename F>
3318 using F_t = std::decay_t<F>;
3319 using ColTypes_t =
typename TTraits::CallableTraits<F_t>::arg_types;
3320 using RetType =
typename TTraits::CallableTraits<F_t>::ret_type;
3321 constexpr auto nColumns = ColTypes_t::list_size;
3336 auto variation = std::make_shared<RDFInternal::RVariation<F_t, IsSingleColumn>>(
3367 throw std::logic_error(
"A column name was passed to the same Vary invocation multiple times.");
3383 template <
typename Helper,
typename ActionResultType>
3385 const std::shared_ptr<Helper> &
hPtr,
3395 const std::shared_ptr<Helper>& ,
3398 throw std::logic_error(std::string(
"An action was booked with no input columns, but the action requires "
3399 "columns! The action helper type was ") +
3400 typeid(Helper).
name());
unsigned long long ULong64_t
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Base class for action helpers, see RInterface::Book() for more information.
implementation of FilterAvailable and FilterMissing operations
The head node of a RDF computation graph.
Helper class that provides the operation graph nodes.
A RDataFrame node that produces a result.
A binder for user-defined columns, variations and aliases.
std::vector< std::string_view > GenerateColumnNames() const
Return the list of the names of the defined columns (Defines + Aliases).
The dataset specification for RDataFrame.
virtual const std::vector< std::string > & GetColumnNames() const =0
Returns a reference to the collection of the dataset's column names.
ColumnNames_t GetValidatedColumnNames(const unsigned int nColumns, const ColumnNames_t &columns)
ColumnNames_t GetColumnTypeNamesList(const ColumnNames_t &columnList)
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > fLoopManager
< The RLoopManager at the root of this computation graph. Never null.
RResultPtr< ActionResultType > CreateAction(const ColumnNames_t &columns, const std::shared_ptr< ActionResultType > &r, const std::shared_ptr< HelperArgType > &helperArg, const std::shared_ptr< RDFNode > &proxiedPtr, const int=-1)
Create RAction object, return RResultPtr for the action Overload for the case in which all column typ...
RDataSource * GetDataSource() const
void CheckAndFillDSColumns(ColumnNames_t validCols, TTraits::TypeList< ColumnTypes... > typeList)
void CheckIMTDisabled(std::string_view callerName)
ColumnNames_t GetColumnNames()
Returns the names of the available columns.
RDFInternal::RColumnRegister fColRegister
Contains the columns defined up to this node.
The public interface to the RDataFrame federation of classes.
RResultPtr<::THnD > HistoND(const THnDModel &model, const ColumnNames_t &columnList)
Fill and return an N-dimensional histogram (lazy action).
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).
RInterface(const std::shared_ptr< Proxied > &proxied, RLoopManager &lm, const RDFInternal::RColumnRegister &colRegister)
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).
RResultPtr<::THnD > HistoND(const THnDModel &model, const ColumnNames_t &columnList)
Fill and return an N-dimensional histogram (lazy action).
std::enable_if_t<!IsFStringConv &&!IsRetTypeDefConstr, RInterface< Proxied, DS_t > > DefineImpl(std::string_view, F, const ColumnNames_t &, const std::string &)
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).
RInterface< Proxied, DS_t > Vary(std::string_view colName, F &&expression, const ColumnNames_t &inputColumns, std::size_t nVariations, std::string_view variationName="")
Register systematic variations for a single existing column using auto-generated variation tags.
RInterface< Proxied, DS_t > Vary(std::string_view colName, std::string_view expression, std::size_t nVariations, std::string_view variationName="")
Register systematic variations for a single existing column using auto-generated variation tags.
RResultPtr<::TGraph > Graph(std::string_view x="", std::string_view y="")
Fill and return a TGraph object (lazy action).
RResultPtr< ActionResultType > CallCreateActionWithoutColsIfPossible(const std::shared_ptr< ActionResultType > &, const std::shared_ptr< Helper > &, Others...)
RInterface< Proxied, DS_t > DefineSlot(std::string_view name, F expression, const ColumnNames_t &columns={})
Define a new 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).
std::enable_if_t< std::is_default_constructible< RetType >::value, RInterface< Proxied, DS_t > > DefineImpl(std::string_view name, F &&expression, const ColumnNames_t &columns, const std::string &where)
RInterface< Proxied, DS_t > DefinePerSample(std::string_view name, F expression)
Define a new column that is updated when the input sample changes.
RInterface & operator=(RInterface &&)=default
Move-assignment operator for RInterface.
RInterface< Proxied, DS_t > Vary(const std::vector< std::string > &colNames, F &&expression, const ColumnNames_t &inputColumns, std::size_t nVariations, std::string_view variationName)
Register systematic variations for multiple existing columns using auto-generated tags.
void ForeachSlot(F f, const ColumnNames_t &columns={})
Execute a user-defined function requiring a processing slot index on each entry (instant action).
RInterface< Proxied, DS_t > Vary(std::string_view colName, std::string_view expression, const std::vector< std::string > &variationTags, std::string_view variationName="")
Register systematic variations for a single existing column using custom variation tags.
RResultPtr< RDisplay > Display(const ColumnNames_t &columnList, size_t nRows=5, size_t nMaxCollectionElements=10)
Provides a representation of the columns in the dataset.
RInterface< RLoopManager > Cache(const ColumnNames_t &columnList)
Save selected columns in memory.
RInterface< Proxied, DS_t > Define(std::string_view name, F expression, const ColumnNames_t &columns={})
Define a new column.
RResultPtr< TStatistic > Stats(std::string_view value, std::string_view weight)
Return a TStatistic object, filled once per event (lazy action).
RInterface< Proxied, DS_t > Redefine(std::string_view name, std::string_view expression)
Overwrite the value and/or type of an existing column.
auto CallCreateActionWithoutColsIfPossible(const std::shared_ptr< ActionResultType > &resPtr, const std::shared_ptr< Helper > &hPtr, TTraits::TypeList< RDFDetail::RInferredType >) -> decltype(hPtr->Exec(0u), RResultPtr< ActionResultType >{})
RInterface< Proxied, DS_t > Vary(const std::vector< std::string > &colNames, std::string_view expression, std::size_t nVariations, std::string_view variationName)
Register systematic variations for multiple existing columns using auto-generated variation tags.
RResultPtr<::TH2D > Histo2D(const TH2DModel &model, std::string_view v1Name="", std::string_view v2Name="")
Fill and return a two-dimensional histogram (lazy action).
RResultPtr< RInterface< RLoopManager > > SnapshotImpl(std::string_view fullTreeName, std::string_view filename, const ColumnNames_t &columnList, const RSnapshotOptions &options)
RInterface< Proxied, DS_t > Vary(std::initializer_list< std::string > colNames, F &&expression, const ColumnNames_t &inputColumns, const std::vector< std::string > &variationTags, std::string_view variationName)
Register systematic variations for multiple existing columns using custom variation tags.
RResultPtr<::TProfile > Profile1D(const TProfile1DModel &model)
Fill and return a one-dimensional profile (lazy action).
RInterface(const std::shared_ptr< RLoopManager > &proxied)
Build a RInterface from a RLoopManager.
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, const std::initializer_list< std::string > &columns)
Append a filter to the call graph.
RInterface< Proxied, DS_t > DefinePerSample(std::string_view name, std::string_view expression)
Define a new column that is updated when the input sample changes.
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.
RResultPtr< RDisplay > Display(std::initializer_list< std::string > columnList, size_t nRows=5, size_t nMaxCollectionElements=10)
Provides a representation of the columns in the dataset.
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< Proxied, DS_t > Redefine(std::string_view name, F expression, const ColumnNames_t &columns={})
Overwrite the value and/or type of an existing column.
RInterface< RLoopManager > Cache(std::string_view columnNameRegexp="")
Save selected columns in memory.
RInterface< Proxied, DS_t > VaryImpl(const std::vector< std::string > &colNames, F &&expression, const ColumnNames_t &inputColumns, const std::vector< std::string > &variationTags, std::string_view variationName)
RResultPtr< typename std::decay_t< 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(std::string_view columnNameRegexp="", size_t nRows=5, size_t nMaxCollectionElements=10)
Provides a representation of the columns in the dataset.
RInterface< RDFDetail::RFilterWithMissingValues< Proxied >, DS_t > FilterAvailable(std::string_view column)
Discard entries with missing values.
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).
RInterface< Proxied, DS_t > Vary(std::string_view colName, F &&expression, const ColumnNames_t &inputColumns, const std::vector< std::string > &variationTags, std::string_view variationName="")
Register systematic variations for a single existing column using custom variation tags.
RResultPtr< ULong64_t > Count()
Return the number of entries processed (lazy action).
RInterface< Proxied, DS_t > Vary(const std::vector< std::string > &colNames, std::string_view expression, const std::vector< std::string > &variationTags, std::string_view variationName)
Register systematic variations for multiple existing columns using custom variation tags.
RInterface< Proxied, DS_t > Define(std::string_view name, std::string_view expression)
Define a new 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).
RInterface< Proxied, DS_t > Vary(const std::vector< std::string > &colNames, F &&expression, const ColumnNames_t &inputColumns, const std::vector< std::string > &variationTags, std::string_view variationName)
Register systematic variations for multiple existing columns using custom variation tags.
RInterface< Proxied, DS_t > RedefineSlotEntry(std::string_view name, F expression, const ColumnNames_t &columns={})
Overwrite the value and/or type of an existing column.
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< RLoopManager > CacheImpl(const ColumnNames_t &columnList, std::index_sequence< S... >)
Implementation of cache.
RInterface< RDFDetail::RRange< Proxied >, DS_t > Range(unsigned int end)
Creates a node that filters entries based on range.
RInterface< RDFDetail::RFilterWithMissingValues< Proxied >, DS_t > FilterMissing(std::string_view column)
Keep only the entries that have missing values.
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).
const std::shared_ptr< Proxied > & GetProxiedPtr() const
RInterface< Proxied, DS_t > JittedVaryImpl(const std::vector< std::string > &colNames, std::string_view expression, const std::vector< std::string > &variationTags, std::string_view variationName, bool isSingleColumn)
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).
RInterface< Proxied, DS_t > Vary(std::initializer_list< std::string > colNames, F &&expression, const ColumnNames_t &inputColumns, std::size_t nVariations, std::string_view variationName)
Register systematic variations for for multiple existing columns using custom variation tags.
RResultPtr< std::decay_t< 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).
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< RDisplay > Display(const ColumnNames_t &columnList, size_t nRows=5, size_t nMaxCollectionElements=10)
Provides a representation of the columns in the dataset.
RResultPtr< RCutFlowReport > Report()
Gather filtering statistics.
RInterface< Proxied, DS_t > RedefineSlot(std::string_view name, F expression, const ColumnNames_t &columns={})
Overwrite the value and/or type of an existing column.
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<::TGraphAsymmErrors > GraphAsymmErrors(std::string_view x="", std::string_view y="", std::string_view exl="", std::string_view exh="", std::string_view eyl="", std::string_view eyh="")
Fill and return a TGraphAsymmErrors object (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={})
Define a new 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)
Fill and return a two-dimensional profile (lazy action).
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.
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).
RInterface< Proxied, DS_t > DefaultValueFor(std::string_view column, const T &defaultValue)
In case the value in the given column is missing, provide a default value.
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 > 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).
RInterface< Proxied, DS_t > Vary(std::initializer_list< std::string > colNames, std::string_view expression, std::size_t nVariations, std::string_view variationName)
Register systematic variations for multiple existing columns using auto-generated variation tags.
A RDataSource implementation which is built on top of result proxies.
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
const_iterator begin() const
const_iterator end() const
typename RemoveFirstParameter< T >::type RemoveFirstParameter_t
TDirectory::TContext keeps track and restore the current directory.
A TGraph is an object made of two arrays X and Y with npoints each.
Statistical variable, defined by its mean and variance (RMS).
void CheckForNoVariations(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister)
Throw if the column has systematic variations attached.
ParsedTreePath ParseTreePath(std::string_view fullTreeName)
void CheckForRedefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is already there.
void CheckForDefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is not already there.
void ChangeEmptyEntryRange(const ROOT::RDF::RNode &node, std::pair< ULong64_t, ULong64_t > &&newRange)
std::shared_ptr< RJittedDefine > BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm, RDataSource *ds, const RColumnRegister &colRegister, const ColumnNames_t &branches, std::shared_ptr< RNodeBase > *upcastNodeOnHeap)
Book the jitting of a Define call.
void CheckValidCppVarName(std::string_view var, const std::string &where)
void ChangeSpec(const ROOT::RDF::RNode &node, ROOT::RDF::Experimental::RDatasetSpec &&spec)
Changes the input dataset specification of an RDataFrame.
const std::vector< std::string > & GetTopLevelFieldNames(const ROOT::RDF::RDataSource &ds)
void RemoveDuplicates(ColumnNames_t &columnNames)
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)
std::string GetDataSourceLabel(const ROOT::RDF::RNode &node)
std::string PrettyPrintAddr(const void *const addr)
void TriggerRun(ROOT::RDF::RNode node)
Trigger the execution of an RDataFrame computation graph.
void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames)
std::string DemangleTypeIdName(const std::type_info &typeInfo)
bool AtLeastOneEmptyString(const std::vector< std::string_view > strings)
std::pair< std::vector< std::string >, std::vector< std::string > > AddSizeBranches(const std::vector< std::string > &branches, ROOT::RDF::RDataSource *ds, std::vector< std::string > &&colsWithoutAliases, std::vector< std::string > &&colsWithAliases)
Return copies of colsWithoutAliases and colsWithAliases with size branches for variable-sized array b...
void SetTTreeLifeline(ROOT::RDF::RNode &node, std::any lifeline)
std::shared_ptr< RDFDetail::RJittedFilter > BookFilterJit(std::shared_ptr< RDFDetail::RNodeBase > *prevNodeOnHeap, std::string_view name, std::string_view expression, const ColumnNames_t &branches, const RColumnRegister &colRegister, TTree *tree, RDataSource *ds)
Book the jitting of a Filter call.
ColumnNames_t FilterArraySizeColNames(const ColumnNames_t &columnNames, const std::string &action)
Take a list of column names, return that list with entries starting by '#' filtered out.
std::shared_ptr< RJittedVariation > BookVariationJit(const std::vector< std::string > &colNames, std::string_view variationName, const std::vector< std::string > &variationTags, std::string_view expression, RLoopManager &lm, RDataSource *ds, const RColumnRegister &colRegister, const ColumnNames_t &branches, std::shared_ptr< RNodeBase > *upcastNodeOnHeap, bool isSingleColumn)
Book the jitting of a Vary call.
void CheckForDuplicateSnapshotColumns(const ColumnNames_t &cols)
ColumnNames_t ConvertRegexToColumns(const ColumnNames_t &colNames, std::string_view columnNameRegexp, std::string_view callerName)
std::shared_ptr< RJittedDefine > BookDefinePerSampleJit(std::string_view name, std::string_view expression, RLoopManager &lm, const RColumnRegister &colRegister, std::shared_ptr< RNodeBase > *upcastNodeOnHeap)
Book the jitting of a DefinePerSample call.
void ChangeBeginAndEndEntries(const RNode &node, Long64_t begin, Long64_t end)
std::vector< std::string > GetTopLevelBranchNames(TTree &t)
Get all the top-level branches names, including the ones of the friend trees.
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
std::vector< std::string > ColumnNames_t
ROOT type_traits extensions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
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.
void DisableImplicitMT()
Disables the implicit multi-threading in ROOT (see EnableImplicitMT).
type is TypeList if MustRemove is false, otherwise it is a TypeList with the first type removed
A collection of options to steer the creation of the dataset on file.
ESnapshotOutputFormat fOutputFormat
Which data format to write to.
bool fLazy
Do not start the event loop when Snapshot is called.
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 THnD.
std::shared_ptr<::THnD > 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.