11#ifndef ROOT_RDF_RINTERFACEBASE
12#define ROOT_RDF_RINTERFACEBASE
33class RVariationsDescription;
70 template <
typename RetType>
71 void SanityChecksForVary(
const std::vector<std::string> &colNames,
const std::vector<std::string> &variationTags,
72 std::string_view variationName)
74 R__ASSERT(variationTags.size() > 0 &&
"Must have at least one variation.");
75 R__ASSERT(colNames.size() > 0 &&
"Must have at least one varied column.");
76 R__ASSERT(!variationName.empty() &&
"Must provide a variation name.");
78 for (
auto &colName : colNames) {
86 if (colNames.size() > 1) {
89 throw std::runtime_error(
"This Vary call is varying multiple columns simultaneously but the expression "
90 "does not return an RVec of RVecs.");
93 auto allColTypesEqual =
94 std::all_of(colTypes.begin() + 1, colTypes.end(), [&](
const std::string &t) { return t == colTypes[0]; });
95 if (!allColTypesEqual)
96 throw std::runtime_error(
"Cannot simultaneously vary multiple columns of different types.");
98 const auto &innerTypeID =
typeid(RDFInternal::InnerValueType_t<RetType>);
100 for (
auto i = 0u; i < colTypes.size(); ++i) {
103 if (innerTypeID != *expectedTypeID)
104 throw std::runtime_error(
"Varied values for column \"" + colNames[i] +
"\" have a different type (" +
109 const auto &retTypeID =
typeid(
typename RetType::value_type);
110 const auto &colName = colNames[0];
112 const auto *expectedTypeID =
114 if (retTypeID != *expectedTypeID)
115 throw std::runtime_error(
"Varied values for column \"" + colName +
"\" have a different type (" +
121 if (colNames.size() > 1) {
122 std::set<std::string> uniqueCols(colNames.begin(), colNames.end());
123 if (uniqueCols.size() != colNames.size())
124 throw std::logic_error(
"A column name was passed to the same Vary invocation multiple times.");
135 template <
typename... ColumnTypes>
147 template <
typename ActionTag,
typename... ColTypes,
typename ActionResultType,
typename RDFNode,
148 typename HelperArgType = ActionResultType,
149 std::enable_if_t<!RDFInternal::RNeedJitting<ColTypes...>
::value,
int> = 0>
151 const std::shared_ptr<HelperArgType> &helperArg,
152 const std::shared_ptr<RDFNode> &proxiedPtr,
const int = -1)
154 constexpr auto nColumns =
sizeof...(ColTypes);
161 auto action = RDFInternal::BuildAction<ColTypes...>(validColumnNames, helperArg, nSlots, proxiedPtr, ActionTag{},
170 template <
typename ActionTag,
typename... ColTypes,
typename ActionResultType,
typename RDFNode,
171 typename HelperArgType = ActionResultType,
172 std::enable_if_t<RDFInternal::RNeedJitting<ColTypes...>
::value,
int> = 0>
174 const std::shared_ptr<HelperArgType> &helperArg,
175 const std::shared_ptr<RDFNode> &proxiedPtr,
const int nColumns = -1)
177 auto realNColumns = (nColumns > -1 ? nColumns :
sizeof...(ColTypes));
183 auto *helperArgOnHeap = RDFInternal::MakeSharedOnHeap(helperArg);
187 const auto jittedAction = std::make_shared<RDFInternal::RJittedAction>(*
fLoopManager, validColumnNames,
189 auto jittedActionOnHeap = RDFInternal::MakeWeakOnHeap(jittedAction);
195 return MakeResultPtr(
r, *
fLoopManager, std::move(jittedAction));
209 bool HasColumn(std::string_view columnName);
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 r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
virtual const std::type_info & GetTypeId() const =0
The head node of a RDF computation graph.
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
void ToJitExec(const std::string &) const
unsigned int GetNSlots() const
A binder for user-defined columns, variations and aliases.
RDFDetail::RDefineBase * GetDefine(const std::string &colName) const
Return the RDefine for the requested column name, or nullptr.
A DFDescription contains useful information about a given RDataFrame computation graph.
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.
RVariationsDescription GetVariations() const
Return a descriptor for the systematic variations registered in this branch of the computation graph.
std::string GetColumnType(std::string_view column)
Return the type of a given column as a string.
ColumnNames_t GetValidatedColumnNames(const unsigned int nColumns, const ColumnNames_t &columns)
RDFDescription Describe()
Return information about the dataframe.
ColumnNames_t GetColumnTypeNamesList(const ColumnNames_t &columnList)
RDFDetail::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...
unsigned int GetNRuns() const
Gets the number of event loops run.
RDataSource * fDataSource
Non-owning pointer to a data-source object. Null if no data-source. RLoopManager has ownership of the...
void SanityChecksForVary(const std::vector< std::string > &colNames, const std::vector< std::string > &variationTags, std::string_view variationName)
void CheckAndFillDSColumns(ColumnNames_t validCols, TTraits::TypeList< ColumnTypes... > typeList)
ColumnNames_t GetDefinedColumnNames()
Returns the names of the defined columns.
void CheckIMTDisabled(std::string_view callerName)
unsigned int GetNSlots() const
Gets the number of data processing slots.
bool HasColumn(std::string_view columnName)
Checks if a column is present in the dataset.
std::string DescribeDataset() const
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 nColumns=-1)
Create RAction object, return RResultPtr for the action Overload for the case in which one or more co...
ColumnNames_t GetColumnNames()
Returns the names of the available columns.
RDFDetail::RLoopManager * GetLoopManager() const
RDFInternal::RColumnRegister fColRegister
Contains the columns defined up to this node.
Smart pointer for the return type of actions.
A descriptor for the systematic variations known to a given RDataFrame node.
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.
const std::type_info & TypeName2TypeID(const std::string &name)
Return the type_info associated to a name.
void CheckValidCppVarName(std::string_view var, const std::string &where)
ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns, const RColumnRegister &colRegister, RDataSource *ds)
Given the desired number of columns and the user-provided list of columns:
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::string JitBuildAction(const ColumnNames_t &cols, std::shared_ptr< RDFDetail::RNodeBase > *prevNode, const std::type_info &helperArgType, const std::type_info &at, void *helperArgOnHeap, TTree *tree, const unsigned int nSlots, const RColumnRegister &colRegister, RDataSource *ds, std::weak_ptr< RJittedAction > *jittedActionOnHeap)
std::vector< std::string > ColumnNames_t
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Lightweight storage for a collection of types.