Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RActionSnapshot.hxx
Go to the documentation of this file.
1// Author: Vincenzo Eduardo Padulano CERN 06/2025
2
3/*************************************************************************
4 * Copyright (C) 1995-2025, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11#ifndef ROOT_RACTIONSNAPSHOT
12#define ROOT_RACTIONSNAPSHOT
13
18
19#include <cstddef> // std::size_t
20#include <memory>
21#include <string>
22#include <vector>
23
24namespace ROOT::Internal::RDF {
25
26namespace GraphDrawing {
27std::shared_ptr<GraphNode> AddDefinesToGraph(std::shared_ptr<GraphNode> node, const RColumnRegister &colRegister,
28 const std::vector<std::string> &prevNodeDefines,
29 std::unordered_map<void *, std::shared_ptr<GraphNode>> &visitedMap);
30} // namespace GraphDrawing
31
32template <typename Helper, typename PrevNode>
34
35 // Template needed to avoid dependency on ActionHelpers.hxx
36 Helper fHelper;
37
38 /// Pointer to the previous node in this branch of the computation graph
39 std::vector<std::shared_ptr<PrevNode>> fPrevNodes;
40
41 /// Column readers per slot and per input column
42 std::vector<std::vector<RColumnReaderBase *>> fValues;
43
44 /// The nth flag signals whether the nth input column is a custom column or not.
45 std::vector<bool> fIsDefine;
46
47 /// Types of the columns to Snapshot
48 std::vector<const std::type_info *> fColTypeIDs;
49
50 ROOT::RDF::SampleCallback_t GetSampleCallback() final { return fHelper.GetSampleCallback(); }
51
52public:
53 RActionSnapshot(Helper &&h, const std::vector<std::string> &columns,
54 const std::vector<const std::type_info *> &colTypeIDs, std::shared_ptr<PrevNode> pd,
56 : RActionBase(pd->GetLoopManagerUnchecked(), columns, colRegister, pd->GetVariations()),
57 fHelper(std::move(h)),
58 fPrevNodes{std::move(pd)},
59 fValues(GetNSlots()),
60 fColTypeIDs(colTypeIDs)
61 {
62 fLoopManager->Register(this);
63
64 const auto nColumns = columns.size();
65 fIsDefine.reserve(nColumns);
66 for (auto i = 0u; i < nColumns; ++i)
67 fIsDefine.push_back(colRegister.IsDefineOrAlias(columns[i]));
68 }
69
74
75 ~RActionSnapshot() final { fLoopManager->Deregister(this); }
76
77 /**
78 Retrieve a wrapper to the result of the action that knows how to merge
79 with others of the same type.
80 */
81 std::unique_ptr<ROOT::Detail::RDF::RMergeableValueBase> GetMergeableValue() const final
82 {
83 return fHelper.GetMergeableValue();
84 }
85
86 void Initialize() final { fHelper.Initialize(); }
87
88 void InitSlot(TTreeReader *r, unsigned int slot) final
89 {
90 fValues[slot] = GetUntypedColumnReaders(slot, r, RActionBase::GetColRegister(), *fLoopManager,
91 RActionBase::GetColumnNames(), fColTypeIDs);
92 fHelper.InitTask(r, slot);
93 }
94
95 void *GetValue(unsigned int slot, std::size_t readerIdx, Long64_t entry)
96 {
97 assert(slot < fValues.size());
98 assert(readerIdx < fValues[slot].size());
99 if (auto *val = fValues[slot][readerIdx]->template TryGet<void>(entry))
100 return val;
101
102 throw std::out_of_range{"RDataFrame: Action (" + fHelper.GetActionName() +
103 ") could not retrieve value for column '" + fColumnNames[readerIdx] + "' for entry " +
104 std::to_string(entry) +
105 ". You can use the DefaultValueFor operation to provide a default value, or "
106 "FilterAvailable/FilterMissing to discard/keep entries with missing values instead."};
107 }
108
109 void CallExec(unsigned int slot, Long64_t entry)
110 {
111 std::vector<void *> untypedValues;
112 auto nReaders = fValues[slot].size();
113 untypedValues.reserve(nReaders);
114 for (decltype(nReaders) readerIdx{}; readerIdx < nReaders; readerIdx++)
115 untypedValues.push_back(GetValue(slot, readerIdx, entry));
116
117 fHelper.Exec(slot, untypedValues);
118 }
119
120 void Run(unsigned int slot, Long64_t entry) final
121 {
122 // check if entry passes all filters
123 if (fPrevNodes.front()->CheckFilters(slot, entry))
124 CallExec(slot, entry);
125 }
126
128 {
129 for (auto const &node : fPrevNodes)
130 node->IncrChildrenCount();
131 }
132
133 /// Clean-up operations to be performed at the end of a task.
134 void FinalizeSlot(unsigned int slot) final
135 {
136 fValues[slot].clear();
137 fHelper.CallFinalizeTask(slot);
138 }
139
140 /// Clean-up and finalize the action result (e.g. merging slot-local results).
141 /// It invokes the helper's Finalize method.
143 {
144 fHelper.Finalize();
145 SetHasRun();
146 }
147
148 std::shared_ptr<GraphDrawing::GraphNode>
149 GetGraph(std::unordered_map<void *, std::shared_ptr<GraphDrawing::GraphNode>> &visitedMap) final
150 {
151 // Action nodes do not need to go through CreateFilterNode: they are never common nodes between multiple branches
152 const auto nodeType = HasRun() ? GraphDrawing::ENodeType::kUsedAction : GraphDrawing::ENodeType::kAction;
153 auto thisNode = std::make_shared<GraphDrawing::GraphNode>(fHelper.GetActionName(), visitedMap.size(), nodeType);
154 visitedMap[(void *)this] = thisNode;
155
156 for (auto const &node : fPrevNodes) {
157 auto prevNode = node->GetGraph(visitedMap);
158 const auto &prevColumns = prevNode->GetDefinedColumns();
159 auto upmostNode = AddDefinesToGraph(thisNode, GetColRegister(), prevColumns, visitedMap);
160
161 thisNode->AddDefinedColumns(GetColRegister().GenerateColumnNames());
162 upmostNode->SetPrevNode(prevNode);
163 }
164 return thisNode;
165 }
166
167 /// This method is invoked to update a partial result during the event loop, right before passing the result to a
168 /// user-defined callback registered via RResultPtr::RegisterCallback
169 void *PartialUpdate(unsigned int slot) final { return fHelper.CallPartialUpdate(slot); }
170
171 [[maybe_unused]] std::unique_ptr<RActionBase> MakeVariedAction(std::vector<void *> && /*results*/) final
172 {
173 // TODO: Probably we also need an untyped RVariedAction
174 throw std::runtime_error("RDataFrame::Snapshot: Snapshot with systematic variations is not supported yet.");
175 }
176
177 /**
178 * \brief Returns a new action with a cloned helper.
179 *
180 * \param[in] newResult The result to be filled by the new action (needed to clone the helper).
181 * \return A unique pointer to the new action.
182 */
183 std::unique_ptr<RActionBase> CloneAction(void *newResult) final
184 {
185 return std::make_unique<RActionSnapshot>(fHelper.CallMakeNew(newResult), GetColumnNames(), fColTypeIDs,
186 fPrevNodes.front(), GetColRegister());
187 }
188};
189
190} // namespace ROOT::Internal::RDF
191
192#endif // ROOT_RACTIONSNAPSHOT
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
RActionSnapshot & operator=(const RActionSnapshot &)=delete
void CallExec(unsigned int slot, Long64_t entry)
ROOT::RDF::SampleCallback_t GetSampleCallback() final
std::vector< std::shared_ptr< PrevNode > > fPrevNodes
Pointer to the previous node in this branch of the computation graph.
std::unique_ptr< RActionBase > MakeVariedAction(std::vector< void * > &&) final
std::vector< bool > fIsDefine
The nth flag signals whether the nth input column is a custom column or not.
std::vector< std::vector< RColumnReaderBase * > > fValues
Column readers per slot and per input column.
void InitSlot(TTreeReader *r, unsigned int slot) final
void * PartialUpdate(unsigned int slot) final
This method is invoked to update a partial result during the event loop, right before passing the res...
void FinalizeSlot(unsigned int slot) final
Clean-up operations to be performed at the end of a task.
std::unique_ptr< RActionBase > CloneAction(void *newResult) final
Returns a new action with a cloned helper.
RActionSnapshot(Helper &&h, const std::vector< std::string > &columns, const std::vector< const std::type_info * > &colTypeIDs, std::shared_ptr< PrevNode > pd, const RColumnRegister &colRegister)
void * GetValue(unsigned int slot, std::size_t readerIdx, Long64_t entry)
RActionSnapshot(RActionSnapshot &&)=delete
std::vector< const std::type_info * > fColTypeIDs
Types of the columns to Snapshot.
void Finalize() final
Clean-up and finalize the action result (e.g.
std::unique_ptr< ROOT::Detail::RDF::RMergeableValueBase > GetMergeableValue() const final
Retrieve a wrapper to the result of the action that knows how to merge with others of the same type.
std::shared_ptr< GraphDrawing::GraphNode > GetGraph(std::unordered_map< void *, std::shared_ptr< GraphDrawing::GraphNode > > &visitedMap) final
RActionSnapshot(const RActionSnapshot &)=delete
RActionSnapshot & operator=(RActionSnapshot &&)=delete
void Run(unsigned int slot, Long64_t entry) final
A binder for user-defined columns, variations and aliases.
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition TTreeReader.h:46
std::shared_ptr< GraphNode > AddDefinesToGraph(std::shared_ptr< GraphNode > node, const RColumnRegister &colRegister, const std::vector< std::string > &prevNodeDefines, std::unordered_map< void *, std::shared_ptr< GraphNode > > &visitedMap)
unsigned int GetNSlots()
Definition RDFUtils.cxx:384
std::vector< RDFDetail::RColumnReaderBase * > GetUntypedColumnReaders(unsigned int slot, TTreeReader *treeReader, ROOT::Internal::RDF::RColumnRegister &colRegister, ROOT::Detail::RDF::RLoopManager &lm, const std::vector< std::string > &colNames, const std::vector< const std::type_info * > &colTypeIDs, const std::string &variationName="nominal")
std::function< void(unsigned int, const ROOT::RDF::RSampleInfo &)> SampleCallback_t
The type of a data-block callback, registered with an RDataFrame computation graph via e....