VTK  9.4.20241221
vtkCellGridResponders.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
54#ifndef vtkCellGridResponders_h
55#define vtkCellGridResponders_h
56
57#include "vtkCellAttributeCalculator.h" // For API.
58#include "vtkCommonDataModelModule.h" // For export macro
59#include "vtkObject.h"
60#include "vtkSmartPointer.h" // For return values
61#include "vtkStringToken.h" // For API.
62#include "vtkTypeName.h" // For RegisterQueryResponder.
63
64#include <unordered_map>
65#include <unordered_set>
66
67VTK_ABI_NAMESPACE_BEGIN
70class vtkCellMetadata;
71
72class VTKCOMMONDATAMODEL_EXPORT vtkCellGridResponders : public vtkObject
73{
74public:
76
79 using TagSet = std::unordered_map<vtkStringToken, std::unordered_set<vtkStringToken>>;
80
84 {
87
88 bool Matches(const TagSet& providedTags) const;
89 };
90
92 void PrintSelf(ostream& os, vtkIndent indent) override;
93
97 template <typename CellType, typename QueryType, typename ResponderType>
98 void RegisterQueryResponder(ResponderType* responder)
99 {
100 vtkStringToken queryTypeKey = vtk::TypeName<QueryType>();
101 vtkStringToken cellTypeKey = vtk::TypeName<CellType>();
102 this->Responders[queryTypeKey][cellTypeKey] = responder;
103 }
104
110 bool Query(vtkCellMetadata* cellType, vtkCellGridQuery* query);
111
120 template <typename CellType, typename CalculatorType>
122 {
123 if (!calculator)
124 {
125 vtkErrorMacro("Could not register null cell-attribute calculator.");
126 return false;
127 }
128
129 vtkStringToken calculatorBaseKey = vtk::TypeName<CalculatorType>();
130 vtkStringToken cellTypeKey = vtk::TypeName<CellType>();
131 if (!dynamic_cast<CalculatorType*>(calculator))
132 {
133 vtkErrorMacro("Could not register cell-attribute calculator "
134 << calculator->GetClassName() << " as it does not inherit " << calculatorBaseKey.Data()
135 << ".");
136 return false;
137 }
138
139 this->Calculators[calculatorBaseKey][cellTypeKey][tags] = calculator;
140 return true;
141 }
142
143 template <typename CellType, typename CalculatorType>
145 {
146 using namespace vtk::literals;
147 if (!calculator)
148 {
149 vtkErrorMacro("Could not register null cell-attribute calculator.");
150 return false;
151 }
152
153 vtkStringToken calculatorBaseKey = vtk::TypeName<CalculatorType>();
154 vtkStringToken cellTypeKey = vtk::TypeName<CellType>();
155 if (!dynamic_cast<CalculatorType*>(calculator))
156 {
157 vtkErrorMacro("Could not register cell-attribute calculator "
158 << calculator->GetClassName() << " as it does not inherit " << calculatorBaseKey.Data()
159 << ".");
160 return false;
161 }
162
163 TagSet tagsIncludingType = tags;
164 tagsIncludingType["Type"_token] = { cellTypeKey };
165 this->CalculatorRegistry[calculatorBaseKey].push_back(
166 CalculatorForTagSet{ tagsIncludingType, calculator });
167 return true;
168 }
169
172 vtkCellMetadata* cellType, vtkCellAttribute* attrib, const TagSet& tags) const;
173
174 template <typename CalculatorType>
176 vtkCellMetadata* cellType, vtkCellAttribute* cellAttribute) const
177 {
178 auto calculatorKey = vtk::TypeName<CalculatorType>();
179 auto baseResult = this->AttributeCalculator(calculatorKey, cellType, cellAttribute);
180 auto result = CalculatorType::SafeDownCast(baseResult);
181 return result;
182 }
183
184 template <typename CalculatorType>
186 vtkCellMetadata* cellType, vtkCellAttribute* attrib, const TagSet& tags) const
187 {
188 auto calculatorKey = vtk::TypeName<CalculatorType>();
189 auto baseResult = this->AttributeCalculator(calculatorKey, cellType, attrib, tags);
190 auto result = CalculatorType::SafeDownCast(baseResult);
191 return result;
192 }
193
202
203 template <typename CacheType>
204 vtkSmartPointer<CacheType> GetCacheDataAs(std::size_t key, bool createIfAbsent = false)
205 {
207 auto rawCache = this->GetCacheData(key);
208 if (!rawCache)
209 {
210 if (createIfAbsent)
211 {
213 this->Caches[key] = castCache;
214 }
215 return castCache;
216 }
217 castCache = dynamic_cast<CacheType*>(rawCache.GetPointer());
218 return castCache;
219 }
220
233 bool SetCacheData(std::size_t key, vtkSmartPointer<vtkObject> value, bool overwrite = false);
234
235protected:
237 ~vtkCellGridResponders() override = default;
238
239 std::unordered_map<vtkStringToken,
240 std::unordered_map<vtkStringToken, vtkSmartPointer<vtkCellGridResponderBase>>>
242
243 std::unordered_map<std::size_t, vtkSmartPointer<vtkObject>> Caches;
244
246 std::unordered_map<vtkStringToken, std::vector<CalculatorForTagSet>> CalculatorRegistry;
247
257 std::unordered_map<vtkStringToken,
258 std::unordered_map<vtkStringToken,
259 std::unordered_map<vtkStringToken, vtkSmartPointer<vtkCellAttributeCalculator>>>>
261
262private:
264 void operator=(const vtkCellGridResponders&) = delete;
265};
266
267VTK_ABI_NAMESPACE_END
268#endif
Perform a per-cell calculation on a vtkCellAttribute.
A function defined over the physical domain of a vtkCellGrid.
Perform an operation on cells in a vtkCellMetadata instance.
Respond to a query on one particular type of cell.
A container that holds objects able to respond to queries specialized for particular vtkCellMetadata ...
std::unordered_map< vtkStringToken, std::unordered_map< vtkStringToken, vtkSmartPointer< vtkCellGridResponderBase > > > Responders
bool Query(vtkCellMetadata *cellType, vtkCellGridQuery *query)
Invoke a responder for the given query and cell type.
std::unordered_map< vtkStringToken, std::vector< CalculatorForTagSet > > CalculatorRegistry
A map from a calculator base class to a set of registered prototypes.
std::unordered_map< std::size_t, vtkSmartPointer< vtkObject > > Caches
vtkSmartPointer< vtkCellAttributeCalculator > AttributeCalculator(vtkStringToken calculatorType, vtkCellMetadata *cellType, vtkCellAttribute *attrib, const TagSet &tags) const
Fetch an instance of an attribute calculator for the given tags.
bool RegisterCalculator(vtkCellAttributeCalculator *calculator, const TagSet &tags)
void RegisterQueryResponder(ResponderType *responder)
Register responder for processing a cell's data.
vtkSmartPointer< CacheType > GetCacheDataAs(std::size_t key, bool createIfAbsent=false)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
std::unordered_map< vtkStringToken, std::unordered_map< vtkStringToken, std::unordered_map< vtkStringToken, vtkSmartPointer< vtkCellAttributeCalculator > > > > Calculators
Nested maps from cell-attribute-query-type to cell-type to cell-attribute-tag-sets to concrete query ...
std::unordered_map< vtkStringToken, std::unordered_set< vtkStringToken > > TagSet
A map of tag names (such as the cell's type-name) to values of the tag accepted or provided for that ...
vtkCellGridResponders()=default
vtkSmartPointer< vtkObject > GetCacheData(std::size_t key)
Return a cache object given a key.
vtkSmartPointer< CalculatorType > AttributeCalculator(vtkCellMetadata *cellType, vtkCellAttribute *cellAttribute) const
~vtkCellGridResponders() override=default
bool RegisterCalculator(vtkStringToken tags, vtkCellAttributeCalculator *calculator)
Register a vtkCellAttributeCalculator subclass.
static vtkCellGridResponders * New()
vtkSmartPointer< CalculatorType > AttributeCalculator(vtkCellMetadata *cellType, vtkCellAttribute *attrib, const TagSet &tags) const
bool SetCacheData(std::size_t key, vtkSmartPointer< vtkObject > value, bool overwrite=false)
Add a cache entry mapping key to value.
Metadata for a particular type of cell (finite element).
a simple class to control print indentation
Definition vtkIndent.h:108
const char * GetClassName() const
Return the class name as a string.
abstract base class for most VTK objects
Definition vtkObject.h:162
Hold a reference to a vtkObjectBase instance.
static vtkSmartPointer< T > New()
Create an instance of a VTK object.
Represent a string by its integer hash.
const std::string & Data() const
Return the string corresponding to the token.
A record of a registered calculator along with the tags it requires and values of those tags it will ...
bool Matches(const TagSet &providedTags) const
vtkSmartPointer< vtkCellAttributeCalculator > CalculatorPrototype