VTK  9.3.1
vtkCellMetadata.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
24 #ifndef vtkCellMetadata_h
25 #define vtkCellMetadata_h
26 
27 #include "vtkCellGridResponders.h" // For ivar.
28 #include "vtkCommonDataModelModule.h" // for export macro
29 #include "vtkNew.h" // for queries
30 #include "vtkObject.h"
31 #include "vtkSmartPointer.h" // for constructor return values
32 #include "vtkStringToken.h" // for vtkStringToken::Hash
33 #include "vtkTypeName.h" // for vtk::TypeName<>()
34 
35 #include <functional>
36 #include <set>
37 #include <unordered_map>
38 
39 VTK_ABI_NAMESPACE_BEGIN
40 class vtkCellGrid;
41 class vtkCellGridQuery;
43 class vtkCellAttribute;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkCellMetadata : public vtkObject
46 {
47 public:
50  using MetadataConstructor = std::function<vtkSmartPointer<vtkCellMetadata>(vtkCellGrid*)>;
51  using ConstructorMap = std::unordered_map<vtkStringToken, MetadataConstructor>;
52 
53  vtkTypeMacro(vtkCellMetadata, vtkObject);
54  void PrintSelf(ostream& os, vtkIndent indent) override;
55 
59  template <typename Subclass>
60  static bool RegisterType()
61  {
62  vtkStringToken name = vtk::TypeName<Subclass>();
63  auto status = vtkCellMetadata::Constructors.insert(std::make_pair(name, [](vtkCellGrid* grid) {
64  auto result = vtkSmartPointer<Subclass>::New();
65  if (result)
66  {
67  result->SetCellGrid(grid);
68  }
69  return result;
70  }));
71  return status.second; // true if insertion occurred.
72  }
73 
74  template <typename Subclass>
76  {
77  vtkStringToken name = vtk::TypeName<Subclass>();
78  auto metadata = vtkCellMetadata::NewInstance(name, grid);
79  vtkSmartPointer<Subclass> result(Subclass::SafeDownCast(metadata));
80  return result;
81  }
82 
83  static vtkSmartPointer<vtkCellMetadata> NewInstance(
84  vtkStringToken className, vtkCellGrid* grid = nullptr);
85 
92  CellTypeId Hash() { return vtkStringToken(this->GetClassName()).GetId(); }
93 
104  virtual bool SetCellGrid(vtkCellGrid* parent);
105 
107  vtkGetObjectMacro(CellGrid, vtkCellGrid);
108 
111  virtual vtkIdType GetNumberOfCells() { return 0; }
112 
120  bool Query(vtkCellGridQuery* query);
121 
129  virtual void ShallowCopy(vtkCellMetadata* vtkNotUsed(other)) {}
130  virtual void DeepCopy(vtkCellMetadata* vtkNotUsed(other)) {}
131 
134 
135 protected:
136  vtkCellMetadata() = default;
137  ~vtkCellMetadata() override;
138 
139  vtkCellGrid* CellGrid{ nullptr };
140 
143 
144 private:
145  vtkCellMetadata(const vtkCellMetadata&) = delete;
146  void operator=(const vtkCellMetadata&) = delete;
147 };
148 
149 VTK_ABI_NAMESPACE_END
150 #endif
Perform an operation on cells in a vtkCellMetadata instance.
vtkStringToken::Hash CellTypeId
abstract base class for most VTK objects
Definition: vtkObject.h:51
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void DeepCopy(vtkCellMetadata *vtkNotUsed(other))
vtkCellMetadata * NewInstance() const
Metadata for a particular type of cell (finite element).
virtual vtkIdType GetNumberOfCells()
Return the number of cells of this type in the parent cell-grid object.
Hold a reference to a vtkObjectBase instance.
Definition: vtkMeta.h:23
int vtkIdType
Definition: vtkType.h:315
static vtkCellGridResponders * GetResponders()
Return the set of registered responder types.
A container that holds objects able to respond to queries specialized for particular vtkCellMetadata ...
static vtkSmartPointer< T > New()
Create an instance of a VTK object.
Hash GetId() const
Return the token's ID (usually its hash but possibly not in the case of collisions).
Visualization data composed of cells of arbitrary type.
Definition: vtkCellGrid.h:41
std::uint32_t Hash
a simple class to control print indentation
Definition: vtkIndent.h:28
const char * GetClassName() const
Return the class name as a string.
A function defined over the physical domain of a vtkCellGrid.
static vtkNew< vtkCellGridResponders > Responders
virtual void ShallowCopy(vtkCellMetadata *vtkNotUsed(other))
Copy other into this instance (which must be of the same type).
represent and manipulate attribute data in a dataset
Represent a string by its integer hash.
static ConstructorMap Constructors
static bool RegisterType()
Register a subclass of vtkCellMetadata.
static vtkSmartPointer< Subclass > NewInstance(vtkCellGrid *grid=nullptr)
std::unordered_map< vtkStringToken, MetadataConstructor > ConstructorMap
CellTypeId Hash()
Return a hash of the cell type.
std::function< vtkSmartPointer< vtkCellMetadata >(vtkCellGrid *)> MetadataConstructor