11#ifndef ROOT_RDF_TINTERFACE
12#define ROOT_RDF_TINTERFACE
50#include <initializer_list>
60#include <unordered_set>
84template <
typename Proxied,
typename DataSource>
116template <
typename Proxied,
typename DataSource =
void>
125 template <
typename T,
typename W>
227 RDFInternal::CheckFilter(
f);
228 using ColTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
229 constexpr auto nColumns = ColTypes_t::list_size;
263 template <
typename F>
348 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
399 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
475 template <
typename F>
506 template <
typename F>
537 constexpr auto where =
"Define";
585 template <
typename F>
604 template <
typename F>
608 "RedefineSlotEntry");
627 constexpr auto where =
"Redefine";
678 template <
typename T>
681 constexpr auto where{
"DefaultValueFor"};
699 auto newColumn = std::make_shared<ROOT::Internal::RDF::RDefaultValueFor<T>>(
868 template <
typename F>
907 template <
typename F>
955 template <
typename F>
977 template <
typename F>
1017 template <
typename F>
1048 template <
typename F>
1230 constexpr auto where =
"Alias";
1246 [[deprecated(
"Snapshot is not any more a template. You can safely remove the template parameters.")]]
1354 bool isRNTuple =
false) ->
const std::type_info * {
1357 }
catch (
const std::runtime_error &err) {
1361 if (std::string(err.what()).find(
"Cannot extract type_info of type") != std::string::npos) {
1364 if (
colTypeName.rfind(
"CLING_UNKNOWN_TYPE", 0) == 0)
1366 std::string
msg{
"No runtime type information is available for column \"" +
colName +
1368 "\". Thus, it cannot be written to disk with Snapshot. Make sure to generate and load "
1369 "ROOT dictionaries for the type of this column."};
1371 throw std::runtime_error(
msg);
1386 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1394 std::vector<const std::type_info *>
colTypeIDs;
1414 "The default Snapshot output data format is TTree, but the input data format is RNTuple. If you "
1415 "want to Snapshot to RNTuple or suppress this warning, set the appropriate fOutputFormat option in "
1416 "RSnapshotOptions. Note that this current default behaviour might change in the future.");
1421 auto newRDF = std::make_shared<RInterface<RLoopManager>>(
1424 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1432 std::vector<const std::type_info *>
colTypeIDs;
1478 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1512 std::initializer_list<std::string>
columnList,
1557 auto staticSeq = std::make_index_sequence<
sizeof...(ColumnTypes)>();
1584 cacheCall <<
"*reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RLoopManager>*>("
1586 <<
") = reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RNodeBase>*>("
1599 cacheCall <<
">(*reinterpret_cast<std::vector<std::string>*>("
1623 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1665 if (
stride == 0 || (end != 0 && end < begin))
1666 throw std::runtime_error(
"Range: stride must be strictly greater than 0 and end must be greater than begin.");
1702 template <
typename F>
1705 using arg_types =
typename TTraits::CallableTraits<
decltype(
f)>::arg_types_nodecay;
1706 using ret_type =
typename TTraits::CallableTraits<
decltype(
f)>::ret_type;
1732 template <
typename F>
1736 constexpr auto nColumns = ColTypes_t::list_size;
1741 using Helper_t = RDFInternal::ForeachSlotHelper<F>;
1783 std::is_default_constructible<T>::value,
1784 "reduce object cannot be default-constructed. Please provide an initialisation value (redIdentity)");
1824 auto cSPtr = std::make_shared<ULong64_t>(0);
1825 using Helper_t = RDFInternal::CountHelper;
1852 template <
typename T,
typename COLL = std::vector<T>>
1860 using Helper_t = RDFInternal::TakeHelper<T, T, COLL>;
1862 auto valuesPtr = std::make_shared<COLL>();
1895 template <
typename V = RDFDetail::RInferredType>
1902 std::shared_ptr<::TH1D>
h(
nullptr);
1904 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
1905 h = model.GetHistogram();
1908 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
1930 template <
typename V = RDFDetail::RInferredType>
1956 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
1963 std::shared_ptr<::TH1D>
h(
nullptr);
1965 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
1969 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
1993 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2013 template <
typename V,
typename W>
2047 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2050 std::shared_ptr<::TH2D>
h(
nullptr);
2052 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2055 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2056 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2093 std::shared_ptr<::TH2D>
h(
nullptr);
2095 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2098 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2099 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2108 template <
typename V1,
typename V2,
typename W>
2146 std::string_view
v3Name =
"")
2148 std::shared_ptr<::TH3D>
h(
nullptr);
2150 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2153 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2154 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2197 std::shared_ptr<::TH3D>
h(
nullptr);
2199 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2202 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2203 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2212 template <
typename V1,
typename V2,
typename V3,
typename W>
2243 std::shared_ptr<::THnD>
h(
nullptr);
2245 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2248 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2250 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2251 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2277 std::shared_ptr<::THnD>
h(
nullptr);
2279 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2282 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2284 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2285 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2317 std::shared_ptr<::THnSparseD>
h(
nullptr);
2319 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2322 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2324 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2325 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2351 std::shared_ptr<::THnSparseD>
h(
nullptr);
2353 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2356 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2358 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2359 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2394 template <
typename X = RDFDetail::RInferredType,
typename Y = RDFDetail::RInferredType>
2397 auto graph = std::make_shared<::TGraph>();
2461 std::string_view
exh =
"", std::string_view
eyl =
"", std::string_view
eyh =
"")
2463 auto graph = std::make_shared<::TGraphAsymmErrors>();
2505 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2509 std::shared_ptr<::TProfile>
h(
nullptr);
2511 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2515 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2516 throw std::runtime_error(
"Profiles with no axes limits are not supported yet.");
2554 std::shared_ptr<::TProfile>
h(
nullptr);
2556 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2560 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2561 throw std::runtime_error(
"Profile histograms with no axes limits are not supported yet.");
2573 template <
typename V1,
typename V2,
typename W>
2609 std::string_view
v2Name =
"", std::string_view
v3Name =
"")
2611 std::shared_ptr<::TProfile2D>
h(
nullptr);
2613 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2617 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2618 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2659 std::shared_ptr<::TProfile2D>
h(
nullptr);
2661 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2665 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2666 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2677 template <
typename V1,
typename V2,
typename V3,
typename W>
2720 auto h = std::make_shared<std::decay_t<T>>(std::forward<T>(model));
2721 if (!RDFInternal::HistoUtils<T>::HasAxisLimits(*
h)) {
2722 throw std::runtime_error(
"The absence of axes limits is not supported yet.");
2743 template <
typename V = RDFDetail::RInferredType>
2747 if (!
value.empty()) {
2751 if (std::is_same<V, RDFDetail::RInferredType>::value) {
2775 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2779 constexpr auto vIsInferred = std::is_same<V, RDFDetail::RInferredType>::value;
2780 constexpr auto wIsInferred = std::is_same<W, RDFDetail::RInferredType>::value;
2789 std::string error(
"The ");
2791 error +=
"column type is explicit, while the ";
2793 error +=
" is specified to be inferred. This case is not supported: please specify both types or none.";
2794 throw std::runtime_error(error);
2821 template <
typename T = RDFDetail::RInferredType>
2825 using RetType_t = RDFDetail::MinReturnType_t<T>;
2826 auto minV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::max());
2851 template <
typename T = RDFDetail::RInferredType>
2855 using RetType_t = RDFDetail::MaxReturnType_t<T>;
2856 auto maxV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::lowest());
2882 template <
typename T = RDFDetail::RInferredType>
2886 auto meanV = std::make_shared<double>(0);
2910 template <
typename T = RDFDetail::RInferredType>
2941 template <
typename T = RDFDetail::RInferredType>
2944 const RDFDetail::SumReturnType_t<T> &
initValue = RDFDetail::SumReturnType_t<T>{})
2947 auto sumV = std::make_shared<RDFDetail::SumReturnType_t<T>>(
initValue);
2987 auto rep = std::make_shared<RCutFlowReport>();
2988 using Helper_t = RDFInternal::ReportHelper<Proxied>;
3059 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
3060 typename ArgTypesNoDecay =
typename TTraits::CallableTraits<AccFun>::arg_types_nodecay,
3061 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
3062 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
3072 using Helper_t = RDFInternal::AggregateHelper<AccFun, MergeFun, R, T, U>;
3074 auto action = std::make_unique<Action_t>(
3094 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
3095 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
3096 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
3100 std::is_default_constructible<U>::value,
3101 "aggregated object cannot be default-constructed. Please provide an initialisation value (aggIdentity)");
3171 using HelperT = std::decay_t<Helper>;
3174 static_assert(std::is_base_of<AH, HelperT>::value && std::is_convertible<HelperT *, AH *>::value,
3175 "Action helper of type T must publicly inherit from ROOT::Detail::RDF::RActionImpl<T>");
3177 auto hPtr = std::make_shared<HelperT>(std::forward<Helper>(
helper));
3180 if (std::is_same<FirstColumn, RDFDetail::RInferredType>::value &&
columns.empty()) {
3287 if (
where.compare(0, 8,
"Redefine") != 0) {
3297 using ArgTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
3299 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::Slot>::value,
ArgTypes_t>
::type;
3301 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::SlotAndEntry>::value,
ColTypesTmp_t>
::type;
3303 constexpr auto nColumns = ColTypes_t::list_size;
3333 bool IsFStringConv = std::is_convertible<F, std::string>::value,
3335 std::enable_if_t<!IsFStringConv && !IsRetTypeDefConstr, RInterface<Proxied, DS_t>>
3338 static_assert(std::is_default_constructible<typename TTraits::CallableTraits<F>::ret_type>
::value,
3339 "Error in `Define`: type returned by expression is not default-constructible");
3345 template <
typename...
ColTypes, std::size_t... S>
3352 RDFInternal::TEvalAnd<std::is_copy_constructible<ColTypes>::value...>
::value;
3353 static_assert(
areCopyConstructible,
"Columns of a type which is not copy constructible cannot be cached yet.");
3366 template <
bool IsSingleColumn,
typename F>
3371 using F_t = std::decay_t<F>;
3372 using ColTypes_t =
typename TTraits::CallableTraits<F_t>::arg_types;
3373 using RetType =
typename TTraits::CallableTraits<F_t>::ret_type;
3374 constexpr auto nColumns = ColTypes_t::list_size;
3389 auto variation = std::make_shared<RDFInternal::RVariation<F_t, IsSingleColumn>>(
3420 throw std::logic_error(
"A column name was passed to the same Vary invocation multiple times.");
3436 template <
typename Helper,
typename ActionResultType>
3438 const std::shared_ptr<Helper> &
hPtr,
3448 const std::shared_ptr<Helper>& ,
3451 throw std::logic_error(std::string(
"An action was booked with no input columns, but the action requires "
3452 "columns! The action helper type was ") +
3453 typeid(Helper).
name());
Basic types used by ROOT and required by TInterpreter.
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int)
long long Long64_t
Portable signed long integer 8 bytes.
unsigned long long ULong64_t
Portable unsigned long integer 8 bytes.
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).
RDFDetail::RDefineBase * GetDefine(std::string_view colName) const
Return the RDefine for the requested column name, or nullptr.
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.
RDFDetail::RLoopManager * GetLoopManager() const
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 or RNTuple 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<::THnSparseD > HistoNSparseD(const THnSparseDModel &model, const ColumnNames_t &columnList)
Fill and return a sparse N-dimensional histogram (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).
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 or RNTuple 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).
RResultPtr<::THnSparseD > HistoNSparseD(const THnSparseDModel &model, const ColumnNames_t &columnList)
Fill and return a sparse N-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())
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 or RNTuple 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)
const std::type_info & TypeName2TypeID(const std::string &name)
Return the type_info associated to a name.
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, 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...
void CheckSnapshotOptionsFormatCompatibility(const ROOT::RDF::RSnapshotOptions &opts)
std::shared_ptr< RDFDetail::RJittedFilter > BookFilterJit(std::shared_ptr< RDFDetail::RNodeBase > *prevNodeOnHeap, std::string_view name, std::string_view expression, const RColumnRegister &colRegister, TTree *tree, RDataSource *ds)
Book the jitting of a Filter call.
void CheckForDefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is not already there.
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(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...
std::string ColumnName2ColumnTypeName(const std::string &colName, TTree *, RDataSource *, RDefineBase *, bool vector2RVec=true)
Return a string containing the type of the given branch.
void RemoveRNTupleSubFields(ColumnNames_t &columnNames)
void SetTTreeLifeline(ROOT::RDF::RNode &node, std::any lifeline)
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, 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 CheckForRedefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is already there.
void ChangeBeginAndEndEntries(const RNode &node, Long64_t begin, Long64_t end)
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
std::vector< std::string > ColumnNames_t
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.
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
Tag to let data sources use the native data type when creating a column reader.
A collection of options to steer the creation of the dataset on disk through Snapshot().
A struct which stores some basic parameters of a TH1D.
std::shared_ptr<::TH1D > GetHistogram() const
A struct which stores some basic parameters of a TH2D.
std::shared_ptr<::TH2D > GetHistogram() const
A struct which stores some basic parameters of a TH3D.
std::shared_ptr<::TH3D > GetHistogram() const
A struct which stores some basic parameters of a THnD.
std::shared_ptr<::THnD > GetHistogram() const
A struct which stores some basic parameters of a THnSparseD.
std::shared_ptr<::THnSparseD > GetHistogram() const
A struct which stores some basic parameters of a TProfile.
std::shared_ptr<::TProfile > GetProfile() const
A struct which stores some basic parameters of a TProfile2D.
std::shared_ptr<::TProfile2D > GetProfile() const
Lightweight storage for a collection of types.