11#ifndef ROOT_RDF_TINTERFACE
12#define ROOT_RDF_TINTERFACE
57#include <initializer_list>
67#include <unordered_set>
91template <
typename Proxied,
typename DataSource>
123template <
typename Proxied,
typename DataSource =
void>
132 template <
typename T,
typename W>
234 RDFInternal::CheckFilter(
f);
235 using ColTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
236 constexpr auto nColumns = ColTypes_t::list_size;
270 template <
typename F>
355 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
406 throw std::runtime_error(
"Unknown column: \"" + std::string(column) +
"\"");
482 template <
typename F>
513 template <
typename F>
544 constexpr auto where =
"Define";
592 template <
typename F>
611 template <
typename F>
615 "RedefineSlotEntry");
634 constexpr auto where =
"Redefine";
685 template <
typename T>
688 constexpr auto where{
"DefaultValueFor"};
706 auto newColumn = std::make_shared<ROOT::Internal::RDF::RDefaultValueFor<T>>(
875 template <
typename F>
914 template <
typename F>
962 template <
typename F>
984 template <
typename F>
1024 template <
typename F>
1055 template <
typename F>
1237 constexpr auto where =
"Alias";
1264 6, 40,
"Snapshot does not need template arguments anymore, you can safely remove them from this function call.")
1352#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 40, 0)
1353 static_assert(
false &&
"Remove information about change of Snapshot defaut compression settings.");
1356 if (
const char *
suppress = std::getenv(
"ROOT_RDF_SNAPSHOT_INFO"))
1357 if (std::strcmp(
suppress,
"0") == 0)
1360 if (std::strcmp(
suppress,
"0") == 0)
1364 <<
"\n\tIn ROOT 6.38, the default compression settings of Snapshot have been changed from 101 (ZLIB with "
1365 "compression level 1, the TTree default) to 505 (ZSTD with compression level 5). This change may result "
1366 "in smaller Snapshot output dataset size by default. In order to suppress this message, set "
1367 "'ROOT_RDF_SNAPSHOT_INFO=0' in your environment or set 'ROOT.RDF.Snapshot.Info: 0' in your .rootrc "
1392 bool isRNTuple =
false) ->
const std::type_info * {
1395 }
catch (
const std::runtime_error &err) {
1399 if (std::string(err.what()).find(
"Cannot extract type_info of type") != std::string::npos) {
1402 if (
colTypeName.rfind(
"CLING_UNKNOWN_TYPE", 0) == 0)
1404 std::string
msg{
"No runtime type information is available for column \"" +
colName +
1406 "\". Thus, it cannot be written to disk with Snapshot. Make sure to generate and load "
1407 "ROOT dictionaries for the type of this column."};
1409 throw std::runtime_error(
msg);
1424 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1432 std::vector<const std::type_info *>
colTypeIDs;
1452 "The default Snapshot output data format is TTree, but the input data format is RNTuple. If you "
1453 "want to Snapshot to RNTuple or suppress this warning, set the appropriate fOutputFormat option in "
1454 "RSnapshotOptions. Note that this current default behaviour might change in the future.");
1459 auto newRDF = std::make_shared<RInterface<RLoopManager>>(
1462 auto snapHelperArgs = std::make_shared<RDFInternal::SnapshotHelperArgs>(RDFInternal::SnapshotHelperArgs{
1470 std::vector<const std::type_info *>
colTypeIDs;
1516 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1550 std::initializer_list<std::string>
columnList,
1595 auto staticSeq = std::make_index_sequence<
sizeof...(ColumnTypes)>();
1622 cacheCall <<
"*reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RLoopManager>*>("
1624 <<
") = reinterpret_cast<ROOT::RDF::RInterface<ROOT::Detail::RDF::RNodeBase>*>("
1637 cacheCall <<
">(*reinterpret_cast<std::vector<std::string>*>("
1661 [](
const std::string &
name) { return name.size() < 13 || name.substr(0, 13) !=
"R_rdf_sizeof_"; });
1703 if (
stride == 0 || (end != 0 && end < begin))
1704 throw std::runtime_error(
"Range: stride must be strictly greater than 0 and end must be greater than begin.");
1740 template <
typename F>
1743 using arg_types =
typename TTraits::CallableTraits<
decltype(
f)>::arg_types_nodecay;
1744 using ret_type =
typename TTraits::CallableTraits<
decltype(
f)>::ret_type;
1770 template <
typename F>
1774 constexpr auto nColumns = ColTypes_t::list_size;
1779 using Helper_t = RDFInternal::ForeachSlotHelper<F>;
1821 std::is_default_constructible<T>::value,
1822 "reduce object cannot be default-constructed. Please provide an initialisation value (redIdentity)");
1862 auto cSPtr = std::make_shared<ULong64_t>(0);
1863 using Helper_t = RDFInternal::CountHelper;
1890 template <
typename T,
typename COLL = std::vector<T>>
1898 using Helper_t = RDFInternal::TakeHelper<T, T, COLL>;
1900 auto valuesPtr = std::make_shared<COLL>();
1933 template <
typename V = RDFDetail::RInferredType>
1940 std::shared_ptr<::TH1D>
h(
nullptr);
1942 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
1943 h = model.GetHistogram();
1946 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
1968 template <
typename V = RDFDetail::RInferredType>
1994 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2001 std::shared_ptr<::TH1D>
h(
nullptr);
2003 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2007 if (
h->GetXaxis()->GetXmax() ==
h->GetXaxis()->GetXmin())
2031 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2051 template <
typename V,
typename W>
2085 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2088 std::shared_ptr<::TH2D>
h(
nullptr);
2090 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2093 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2094 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2131 std::shared_ptr<::TH2D>
h(
nullptr);
2133 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2136 if (!RDFInternal::HistoUtils<::TH2D>::HasAxisLimits(*
h)) {
2137 throw std::runtime_error(
"2D histograms with no axes limits are not supported yet.");
2146 template <
typename V1,
typename V2,
typename W>
2184 std::string_view
v3Name =
"")
2186 std::shared_ptr<::TH3D>
h(
nullptr);
2188 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2191 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2192 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2235 std::shared_ptr<::TH3D>
h(
nullptr);
2237 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2240 if (!RDFInternal::HistoUtils<::TH3D>::HasAxisLimits(*
h)) {
2241 throw std::runtime_error(
"3D histograms with no axes limits are not supported yet.");
2250 template <
typename V1,
typename V2,
typename V3,
typename W>
2281 std::shared_ptr<::THnD>
h(
nullptr);
2283 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2286 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2288 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2289 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2315 std::shared_ptr<::THnD>
h(
nullptr);
2317 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2320 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2322 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2323 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2355 std::shared_ptr<::THnSparseD>
h(
nullptr);
2357 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2360 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2362 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2363 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2389 std::shared_ptr<::THnSparseD>
h(
nullptr);
2391 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2394 if (
int(
columnList.size()) == (
h->GetNdimensions() + 1)) {
2396 }
else if (
int(
columnList.size()) !=
h->GetNdimensions()) {
2397 throw std::runtime_error(
"Wrong number of columns for the specified number of histogram axes.");
2432 template <
typename X = RDFDetail::RInferredType,
typename Y = RDFDetail::RInferredType>
2435 auto graph = std::make_shared<::TGraph>();
2499 std::string_view
exh =
"", std::string_view
eyl =
"", std::string_view
eyh =
"")
2501 auto graph = std::make_shared<::TGraphAsymmErrors>();
2543 template <
typename V1 = RDFDetail::RInferredType,
typename V2 = RDFDetail::RInferredType>
2547 std::shared_ptr<::TProfile>
h(
nullptr);
2549 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2553 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2554 throw std::runtime_error(
"Profiles with no axes limits are not supported yet.");
2592 std::shared_ptr<::TProfile>
h(
nullptr);
2594 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2598 if (!RDFInternal::HistoUtils<::TProfile>::HasAxisLimits(*
h)) {
2599 throw std::runtime_error(
"Profile histograms with no axes limits are not supported yet.");
2611 template <
typename V1,
typename V2,
typename W>
2647 std::string_view
v2Name =
"", std::string_view
v3Name =
"")
2649 std::shared_ptr<::TProfile2D>
h(
nullptr);
2651 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2655 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2656 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2697 std::shared_ptr<::TProfile2D>
h(
nullptr);
2699 ROOT::Internal::RDF::RIgnoreErrorLevelRAII
iel(
kError);
2703 if (!RDFInternal::HistoUtils<::TProfile2D>::HasAxisLimits(*
h)) {
2704 throw std::runtime_error(
"2D profiles with no axes limits are not supported yet.");
2715 template <
typename V1,
typename V2,
typename V3,
typename W>
2758 auto h = std::make_shared<std::decay_t<T>>(std::forward<T>(model));
2759 if (!RDFInternal::HistoUtils<T>::HasAxisLimits(*
h)) {
2760 throw std::runtime_error(
"The absence of axes limits is not supported yet.");
2781 template <
typename V = RDFDetail::RInferredType>
2785 if (!
value.empty()) {
2789 if (std::is_same<V, RDFDetail::RInferredType>::value) {
2813 template <
typename V = RDFDetail::RInferredType,
typename W = RDFDetail::RInferredType>
2817 constexpr auto vIsInferred = std::is_same<V, RDFDetail::RInferredType>::value;
2818 constexpr auto wIsInferred = std::is_same<W, RDFDetail::RInferredType>::value;
2827 std::string error(
"The ");
2829 error +=
"column type is explicit, while the ";
2831 error +=
" is specified to be inferred. This case is not supported: please specify both types or none.";
2832 throw std::runtime_error(error);
2859 template <
typename T = RDFDetail::RInferredType>
2863 using RetType_t = RDFDetail::MinReturnType_t<T>;
2864 auto minV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::max());
2889 template <
typename T = RDFDetail::RInferredType>
2893 using RetType_t = RDFDetail::MaxReturnType_t<T>;
2894 auto maxV = std::make_shared<RetType_t>(std::numeric_limits<RetType_t>::lowest());
2920 template <
typename T = RDFDetail::RInferredType>
2924 auto meanV = std::make_shared<double>(0);
2948 template <
typename T = RDFDetail::RInferredType>
2979 template <
typename T = RDFDetail::RInferredType>
2982 const RDFDetail::SumReturnType_t<T> &
initValue = RDFDetail::SumReturnType_t<T>{})
2985 auto sumV = std::make_shared<RDFDetail::SumReturnType_t<T>>(
initValue);
3025 auto rep = std::make_shared<RCutFlowReport>();
3026 using Helper_t = RDFInternal::ReportHelper<Proxied>;
3097 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
3098 typename ArgTypesNoDecay =
typename TTraits::CallableTraits<AccFun>::arg_types_nodecay,
3099 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
3100 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
3110 using Helper_t = RDFInternal::AggregateHelper<AccFun, MergeFun, R, T, U>;
3112 auto action = std::make_unique<Action_t>(
3132 typename ArgTypes =
typename TTraits::CallableTraits<AccFun>::arg_types,
3133 typename U = TTraits::TakeFirstParameter_t<ArgTypes>,
3134 typename T = TTraits::TakeFirstParameter_t<TTraits::RemoveFirstParameter_t<ArgTypes>>>
3138 std::is_default_constructible<U>::value,
3139 "aggregated object cannot be default-constructed. Please provide an initialisation value (aggIdentity)");
3209 using HelperT = std::decay_t<Helper>;
3212 static_assert(std::is_base_of<AH, HelperT>::value && std::is_convertible<HelperT *, AH *>::value,
3213 "Action helper of type T must publicly inherit from ROOT::Detail::RDF::RActionImpl<T>");
3215 auto hPtr = std::make_shared<HelperT>(std::forward<Helper>(
helper));
3218 if (std::is_same<FirstColumn, RDFDetail::RInferredType>::value &&
columns.empty()) {
3325 if (
where.compare(0, 8,
"Redefine") != 0) {
3335 using ArgTypes_t =
typename TTraits::CallableTraits<F>::arg_types;
3337 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::Slot>::value,
ArgTypes_t>
::type;
3339 std::is_same<DefineType, RDFDetail::ExtraArgsForDefine::SlotAndEntry>::value,
ColTypesTmp_t>
::type;
3341 constexpr auto nColumns = ColTypes_t::list_size;
3371 bool IsFStringConv = std::is_convertible<F, std::string>::value,
3373 std::enable_if_t<!IsFStringConv && !IsRetTypeDefConstr, RInterface<Proxied, DS_t>>
3376 static_assert(std::is_default_constructible<typename TTraits::CallableTraits<F>::ret_type>
::value,
3377 "Error in `Define`: type returned by expression is not default-constructible");
3383 template <
typename...
ColTypes, std::size_t... S>
3390 RDFInternal::TEvalAnd<std::is_copy_constructible<ColTypes>::value...>
::value;
3391 static_assert(
areCopyConstructible,
"Columns of a type which is not copy constructible cannot be cached yet.");
3404 template <
bool IsSingleColumn,
typename F>
3409 using F_t = std::decay_t<F>;
3410 using ColTypes_t =
typename TTraits::CallableTraits<F_t>::arg_types;
3411 using RetType =
typename TTraits::CallableTraits<F_t>::ret_type;
3412 constexpr auto nColumns = ColTypes_t::list_size;
3427 auto variation = std::make_shared<RDFInternal::RVariation<F_t, IsSingleColumn>>(
3458 throw std::logic_error(
"A column name was passed to the same Vary invocation multiple times.");
3474 template <
typename Helper,
typename ActionResultType>
3476 const std::shared_ptr<Helper> &
hPtr,
3486 const std::shared_ptr<Helper>& ,
3489 throw std::logic_error(std::string(
"An action was booked with no input columns, but the action requires "
3490 "columns! The action helper type was ") +
3491 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).
R__DEPRECATED(6, 40, "Snapshot does not need template arguments anymore, you can safely remove them from this function call.") RResultPtr< RInterface< RLoopManager > > Snapshot(std
Save selected columns to disk, in a new TTree or RNTuple 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 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.
Smart pointer for the return type of actions.
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
Change the verbosity level (global or specific to the RLogChannel passed to the constructor) for the ...
const_iterator begin() const
const_iterator end() const
typename RemoveFirstParameter< T >::type RemoveFirstParameter_t
TDirectory::TContext keeps track and restore the current directory.
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
A TGraph is an object made of two arrays X and Y with npoints each.
Statistical variable, defined by its mean and variance (RMS).
ROOT::RLogChannel & RDFLogChannel()
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.
@ kInfo
Informational messages; used for instance for tracing.
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.