40#pragma GCC diagnostic push 
   41#pragma GCC diagnostic ignored "-Woverloaded-virtual" 
   42#pragma GCC diagnostic ignored "-Wshadow" 
   46#pragma GCC diagnostic pop 
   61#include <unordered_map> 
   62#include <unordered_set> 
   85struct ParsedExpression {
 
   90   ColumnNames_t fUsedCols;
 
   92   ColumnNames_t fVarNames;
 
   96std::pair<ColumnNames_t, ColumnNames_t>
 
  100   lexertk::generator 
tokens;
 
  103      const auto msg = 
"Failed to tokenize expression:\n" + 
expr + 
"\n\nMake sure it is valid C++.";
 
  104      throw std::runtime_error(
msg);
 
  107   std::unordered_set<std::string> 
usedCols;
 
  112   const auto kSymbol = lexertk::token::e_symbol;
 
  113   for (
auto i = 0
u; i < 
nTokens; ++i) {
 
  163   dot.Substitute(out, 
"\\.", 
"g");
 
  164   return std::string(std::move(out));
 
  175      replacer.Substitute(out, col.data(), 
"g");
 
  189      "(^|\\W)#(?!(ifdef|ifndef|if|else|elif|endif|pragma|define|undef|include|line))([a-zA-Z_][a-zA-Z0-9_]*)");
 
  204   for (
auto i = 0
u; i < 
varNames.size(); ++i)
 
  205      varNames[i] = 
"var" + std::to_string(i);
 
  213             [](
const std::string &
a, 
const std::string &
b) { return a.size() > b.size(); });
 
  276   std::stringstream 
ss;
 
  278   for (
auto i = 0
u; i < vars.size(); ++i) {
 
  292      ss.seekp(-2, 
ss.cur);
 
  314      const auto funcName = 
exprIt->second;
 
  323                          "_ret_t = typename ROOT::TypeTraits::CallableTraits<decltype(" + 
funcBaseName +
 
  337   const auto dt = 
gROOT->GetType((funcName + 
"_ret_t").c_str());
 
  339   const auto type = 
dt->GetFullTypeName();
 
  354      "RDataFrame::Jit: cannot just-in-time compile a \"" + 
actionTypeNameBase + 
"\" action using helper type \"" +
 
  356      "\". This typically happens in a custom `Fill` or `Book` invocation where the types of the input columns have " 
  357      "not been specified as template parameters and the ROOT interpreter has no knowledge of this type of action " 
  358      "helper. Please add template parameters for the types of the input columns to avoid jitting this action (i.e. " 
  359      "`df.Fill<float>(..., {\"x\"})`, where `float` is the type of `x`) or declare the action helper type to the " 
  360      "interpreter, e.g. via gInterpreter->Declare.";
 
  378                [&](
const std::string &
name) {
 
  379                   if (name[0] == 
'#') {
 
  380                     filteredColumns.emplace_back(name);
 
  388      std::string 
msg = 
"Column name(s) {";
 
  391      msg[
msg.size() - 2] = 
'}';
 
  392      msg += 
"will be ignored. Please go through a valid Alias to " + 
action + 
" an array size column";
 
  393      throw std::runtime_error(
msg);
 
 
  406   if (col.size() > 1 && col[0] == 
'#')
 
  407      return "R_rdf_sizeof_" + col.substr(1);
 
 
  421   auto isALetter = [](
char c) { 
return (
c >= 
'A' && 
c <= 
'Z') || (
c >= 
'a' && 
c <= 
'z'); };
 
  427   auto isANumber = [](
char c) { 
return c >= 
'0' && 
c <= 
'9'; };
 
  429   for (
const char c : var)
 
  434      const auto objName = 
where == 
"Define" ? 
"column" : 
"variation";
 
  435      const auto error = 
"RDataFrame::" + 
where + 
": cannot define " + 
objName + 
" \"" + std::string(var) +
 
  436                         "\". Not a valid C++ variable name.";
 
  437      throw std::runtime_error(error);
 
 
  441std::string DemangleTypeIdName(
const std::type_info &
typeInfo)
 
 
  477         text = 
": there is no column available to match.";
 
  481      throw std::runtime_error(
text);
 
 
  494              "\", already exists in this branch of the computation graph.";
 
  496      error = 
"A column with that name has already been Define'd. Use Redefine to force redefinition.";
 
  502         "A branch with that name is already present in the input TTree/TChain. Use Redefine to force redefinition.";
 
  505         "A column with that name is already present in the input data source. Use Redefine to force redefinition.";
 
  507   if (!error.empty()) {
 
  508      error = 
"RDataFrame::" + 
where + 
": cannot define column \"" + std::string(
definedColView) + 
"\". " + error;
 
  509      throw std::runtime_error(error);
 
 
  521              "\", already exists. Aliases cannot be Redefined or Varied.";
 
  534         error = 
"No column with that name was found in the dataset. Use Define to create a new column.";
 
  537   if (!error.empty()) {
 
  538      if (
where == 
"DefaultValueFor")
 
  539         error = 
"RDataFrame::" + 
where + 
": cannot provide default values for column \"" +
 
  542         error = 
"RDataFrame::" + 
where + 
": cannot redefine or vary column \"" + std::string(
definedColView) + 
"\". " +
 
  544      throw std::runtime_error(error);
 
 
  554      if (
where == 
"Redefine") {
 
  555         const std::string error = 
"RDataFrame::" + 
where + 
": cannot redefine column \"" + 
definedCol +
 
  556                                   "\". The column depends on one or more systematic variations and re-defining varied " 
  557                                   "columns is not supported.";
 
  558         throw std::runtime_error(error);
 
  559      } 
else if (
where == 
"DefaultValueFor") {
 
  560         const std::string error = 
"RDataFrame::" + 
where + 
": cannot provide a default value for column \"" +
 
  562                                   "\". The column depends on one or more systematic variations and it should not be " 
  563                                   "possible to have missing values in varied columns.";
 
  564         throw std::runtime_error(error);
 
  566         const std::string error =
 
  567            "RDataFrame::" + 
where + 
": this operation cannot work with columns that depend on systematic variations.";
 
  568         throw std::runtime_error(error);
 
 
  576      std::string 
err_msg = 
"The number of template parameters specified is ";
 
  580      err_msg += 
" columns have been specified.";
 
  581      throw std::runtime_error(
err_msg);
 
 
  592         throw std::runtime_error(
 
  594            " required but none were provided and the default list has size " + std::to_string(
defaultNames.size()));
 
  601                    " required but " + std::to_string(names.size()) + (names.size() == 1 ? 
" was" : 
" were") +
 
  603         for (
const auto &
name : names)
 
  606         throw std::runtime_error(
msg);
 
 
  642   if (std::string_view::npos != 
lastSlash) {
 
 
  658std::shared_ptr<RDFDetail::RJittedFilter>
 
  670      std::runtime_error(
"Filter: the following expression does not evaluate to bool:\n" + std::string(expression));
 
  677   const auto jittedFilter = std::make_shared<RDFDetail::RJittedFilter>(
 
  678      (*prevNodeOnHeap)->GetLoopManagerUnchecked(), 
name,
 
  684   filterInvocation << 
"ROOT::Internal::RDF::JitFilterHelper(" << funcName << 
", new const char*[" 
  695                    << 
"reinterpret_cast<std::weak_ptr<ROOT::Detail::RDF::RJittedFilter>*>(" 
  697                    << 
"reinterpret_cast<std::shared_ptr<ROOT::Detail::RDF::RNodeBase>*>(" << 
prevNodeAddr << 
")," 
  698                    << 
"reinterpret_cast<ROOT::Internal::RDF::RColumnRegister*>(" << 
definesOnHeapAddr << 
")" 
 
  713   auto *
const tree = 
lm.GetTree();
 
  727   defineInvocation << 
"ROOT::Internal::RDF::JitDefineHelper<ROOT::Internal::RDF::DefineTypes::RDefineTag>(" << funcName
 
  728                    << 
", new const char*[" << 
parsedExpr.fUsedCols.size() << 
"]{";
 
  729   for (
const auto &col : 
parsedExpr.fUsedCols) {
 
  739                    << 
"\", reinterpret_cast<ROOT::Detail::RDF::RLoopManager*>(" << 
PrettyPrintAddr(&
lm)
 
  740                    << 
"), reinterpret_cast<std::weak_ptr<ROOT::Detail::RDF::RJittedDefine>*>(" 
  742                    << 
"), reinterpret_cast<ROOT::Internal::RDF::RColumnRegister*>(" << 
definesAddr 
  743                    << 
"), reinterpret_cast<std::shared_ptr<ROOT::Detail::RDF::RNodeBase>*>(" 
 
  755   const auto funcName = 
DeclareFunction(std::string(expression), {
"rdfslot_", 
"rdfsampleinfo_"},
 
  756                                         {
"unsigned int", 
"const ROOT::RDF::RSampleInfo"});
 
  764   defineInvocation << 
"ROOT::Internal::RDF::JitDefineHelper<ROOT::Internal::RDF::DefineTypes::RDefinePerSampleTag>(" 
  765                    << funcName << 
", nullptr, 0, ";
 
  771                    << 
"), reinterpret_cast<std::weak_ptr<ROOT::Detail::RDF::RJittedDefine>*>(" 
  773                    << 
"), reinterpret_cast<ROOT::Internal::RDF::RColumnRegister*>(" << 
definesAddr 
  774                    << 
"), reinterpret_cast<std::shared_ptr<ROOT::Detail::RDF::RNodeBase>*>(" 
 
  782std::shared_ptr<RJittedVariation>
 
  788   auto *
const tree = 
lm.GetTree();
 
  797   if (
type.rfind(
"ROOT::VecOps::RVec", 0) != 0) {
 
  801      throw std::runtime_error(
 
  802         "Jitted Vary expressions must return an RVec object. The following expression returns a " + 
type +
 
  819                  << funcName << 
", new const char*[" << 
parsedExpr.fUsedCols.size() << 
"]{";
 
  820   for (
const auto &col : 
parsedExpr.fUsedCols) {
 
  837                  << 
"\", reinterpret_cast<ROOT::Detail::RDF::RLoopManager*>(" << 
PrettyPrintAddr(&
lm)
 
  838                  << 
"), reinterpret_cast<std::weak_ptr<ROOT::Internal::RDF::RJittedVariation>*>(" 
  840                  << 
"), reinterpret_cast<ROOT::Internal::RDF::RColumnRegister*>(" << 
colRegisterAddr 
  841                  << 
"), reinterpret_cast<std::shared_ptr<ROOT::Detail::RDF::RNodeBase>*>(" 
 
  858      std::string 
exceptionText = 
"An error occurred while inferring the action type of the operation.";
 
  882   createAction_str << 
">(reinterpret_cast<std::shared_ptr<ROOT::Detail::RDF::RNodeBase>*>(" 
  884   for (
auto i = 0
u; i < 
cols.size(); ++i) {
 
  891                    << 
"), reinterpret_cast<std::weak_ptr<ROOT::Internal::RDF::RJittedAction>*>(" 
  893                    << 
"), reinterpret_cast<ROOT::Internal::RDF::RColumnRegister*>(" << 
definesAddr << 
"));";
 
 
  899   for (
const auto &s : 
strings) {
 
 
  906std::shared_ptr<RNodeBase> 
UpcastNode(std::shared_ptr<RNodeBase> ptr)
 
 
  935      const auto &
colsToIgnore = 
lm.GetSuppressErrorsForMissingBranches();
 
  940         std::string 
errMsg = std::string(
"Unknown column") + (
unknownColumns.size() > 1 ? 
"s: " : 
": ");
 
  944         throw std::runtime_error(
errMsg);
 
 
  958      if (
colType.rfind(
"CLING_UNKNOWN_TYPE", 0) == 0) { 
 
  960            "The type of custom column \"" + 
c + 
"\" (" + 
colType.substr(19) +
 
  961            ") is not known to the interpreter, but a just-in-time-compiled " + context +
 
  962            " call requires this column. Make sure to create and load ROOT dictionaries for this column's class.";
 
  963         throw std::runtime_error(
msg);
 
 
  988   for (
auto &col : 
cols) {
 
  990         const auto msg = 
"Error: column \"" + col +
 
  991                          "\" was passed to Snapshot twice. This is not supported: only one of the columns would be " 
  992                          "readable with RDataFrame.";
 
  993         throw std::logic_error(
msg);
 
 
 1000std::pair<std::vector<std::string>, std::vector<std::string>>
 
 1006      tree = 
treeDS->GetTree();
 
 1014   for (std::size_t i = 0
u; i < 
nCols; ++i) {
 
 1019      auto *
b = tree->GetBranch(
colName.c_str());
 
 1021         b = tree->FindBranch(
colName.c_str());
 
 1023      auto *
leaves = 
b->GetListOfLeaves();
 
 
 1046                     [&
uniqueCols](
const std::string &
colName) { return !uniqueCols.insert(colName).second; }),
 
 
 1055                [](
const std::string &
colName) { return colName.find(
'.') == std::string::npos; });
 
 1059                                       if (colName.find(
'.') == std::string::npos)
 
 1061                                       const auto parentFieldName = colName.substr(0, colName.find_first_of(
'.'));
 
 1062                                       return std::find(parentFields.cbegin(), parentFields.cend(), parentFieldName) !=
 
 
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
 
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.
 
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 cname
 
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
 
Option_t Option_t TPoint TPoint const char text
 
R__EXTERN TVirtualMutex * gROOTMutex
 
#define R__LOCKGUARD(mutex)
 
The head node of a RDF computation graph.
 
A binder for user-defined columns, variations and aliases.
 
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
 
const_iterator begin() const
 
const_iterator end() const
 
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
 
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
 
Bool_t MatchB(const TString &s, const TString &mods="", Int_t start=0, Int_t nMaxMatch=10)
 
A TTree represents a columnar dataset.
 
const ColumnNames_t SelectColumns(unsigned int nRequiredNames, const ColumnNames_t &names, const ColumnNames_t &defaultNames)
Choose between local column names or default column names, throw in case of errors.
 
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.
 
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 RemoveDuplicates(ColumnNames_t &columnNames)
 
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...
 
bool IsStrInVec(const std::string &str, const std::vector< std::string > &vec)
 
std::string ResolveAlias(const std::string &col, const std::map< std::string, std::string > &aliasMap)
 
std::vector< std::string > GetFilterNames(const std::shared_ptr< RLoopManager > &loopManager)
 
std::string PrettyPrintAddr(const void *const addr)
 
void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames)
 
bool AtLeastOneEmptyString(const std::vector< std::string_view > strings)
 
std::string ColumnName2ColumnTypeName(const std::string &colName, TTree *, RDataSource *, RDefineBase *, bool vector2RVec=true)
Return a string containing the type of the given branch.
 
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 RemoveRNTupleSubFields(ColumnNames_t &columnNames)
 
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.
 
std::vector< T > Union(const std::vector< T > &v1, const std::vector< T > &v2)
Return a vector with all elements of v1 and v2 and duplicates removed.
 
bool IsInternalColumn(std::string_view colName)
Whether custom column with name colName is an "internal" column such as rdfentry_ or rdfslot_.
 
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.
 
void InterpreterDeclare(const std::string &code)
Declare code in the interpreter via the TInterpreter::Declare method, throw in case of errors.
 
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.
 
std::vector< std::string > GetValidatedArgTypes(const ColumnNames_t &colNames, const RColumnRegister &colRegister, TTree *tree, RDataSource *ds, const std::string &context, bool vector2RVec)
 
void CheckForDuplicateSnapshotColumns(const ColumnNames_t &cols)
 
ColumnNames_t ConvertRegexToColumns(const ColumnNames_t &colNames, std::string_view columnNameRegexp, std::string_view callerName)
 
std::vector< bool > FindUndefinedDSColumns(const ColumnNames_t &requestedCols, const ColumnNames_t &definedCols)
Return a bitset each element of which indicates whether the corresponding element in selectedColumns ...
 
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.
 
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, const bool vector2RVec)
 
ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns, const RColumnRegister &definedCols, const ColumnNames_t &dataSourceColumns)
 
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...
 
char * DemangleTypeIdName(const std::type_info &ti, int &errorCode)
Demangle in a portable way the type id name.
 
BVH_ALWAYS_INLINE T dot(const Vec< T, N > &a, const Vec< T, N > &b)