VTK  9.4.20241221
vtkMappedUnstructuredGrid.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
131#ifndef vtkMappedUnstructuredGrid_h
132#define vtkMappedUnstructuredGrid_h
133
135
136#include "vtkMappedUnstructuredGridCellIterator.h" // For default cell iterator
137#include "vtkNew.h" // For vtkNew
138#include "vtkSmartPointer.h" // For vtkSmartPointer
139
140VTK_ABI_NAMESPACE_BEGIN
141template <class Implementation,
144{
146
147public:
149 typedef Implementation ImplementationType;
150 typedef CellIterator CellIteratorType;
151
152 // Virtuals from various base classes:
153 void PrintSelf(ostream& os, vtkIndent indent) override;
154 void CopyStructure(vtkDataSet* pd) override;
155 void ShallowCopy(vtkDataObject* src) override;
158 vtkCell* GetCell(vtkIdType cellId) override;
159 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
160 int GetCellType(vtkIdType cellId) override;
162 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
164 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override;
165 int GetMaxCellSize() override;
166 void GetIdsOfCellsOfType(int type, vtkIdTypeArray* array) override;
167 int IsHomogeneous() override;
168 void Allocate(vtkIdType numCells, int extSize = 1000) override;
170
173
174protected:
177
178 // For convenience...
180
182
183 vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[]) override;
184 vtkIdType InternalInsertNextCell(int type, vtkIdList* ptIds) override;
186 int type, vtkIdType npts, const vtkIdType ptIds[], vtkCellArray* faces) override;
187 void InternalReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) override;
188
189private:
191 void operator=(const vtkMappedUnstructuredGrid&) = delete;
192
193 vtkNew<vtkGenericCell> TempCell;
194};
195
196VTK_ABI_NAMESPACE_END
197#include "vtkMappedUnstructuredGrid.txx"
198
199// We need to fake the superclass for the wrappers, otherwise they will choke on
200// the template:
201#ifndef __VTK_WRAP__
202
203#define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
204 class _exportDecl _className : public vtkMappedUnstructuredGrid<_impl> \
205 { \
206 public: \
207 vtkTypeMacro(_className, vtkMappedUnstructuredGrid<_impl>); \
208 static _className* New(); \
209 \
210 protected: \
211 _className() \
212 { \
213 _impl* i = _impl::New(); \
214 this->SetImplementation(i); \
215 i->Delete(); \
216 } \
217 ~_className() override {} \
218 \
219 private: \
220 _className(const _className&); \
221 void operator=(const _className&); \
222 }
223
224#define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
225 class _exportDecl _className : public vtkMappedUnstructuredGrid<_impl, _cIter> \
226 { \
227 public: \
228 typedef vtkMappedUnstructuredGrid<_impl, _cIter> SelfType; \
229 vtkTypeMacro(_className, SelfType); \
230 static _className* New(); \
231 \
232 protected: \
233 _className() \
234 { \
235 _impl* i = _impl::New(); \
236 this->SetImplementation(i); \
237 i->Delete(); \
238 } \
239 ~_className() override {} \
240 \
241 private: \
242 _className(const _className&); \
243 void operator=(const _className&); \
244 }
245
246#else // __VTK_WRAP__
247
248#define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
249 class _exportDecl _className : public vtkUnstructuredGridBase \
250 { \
251 public: \
252 vtkTypeMacro(_className, vtkUnstructuredGridBase); \
253 static _className* New(); \
254 \
255 protected: \
256 _className() {} \
257 ~_className() override {} \
258 \
259 private: \
260 _className(const _className&); \
261 void operator=(const _className&); \
262 }
263
264#define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
265 class _exportDecl _className : public vtkUnstructuredGridBase \
266 { \
267 public: \
268 vtkTypeMacro(_className, vtkUnstructuredGridBase); \
269 static _className* New(); \
270 \
271 protected: \
272 _className() {} \
273 ~_className() override {} \
274 \
275 private: \
276 _className(const _className&); \
277 void operator=(const _className&); \
278 }
279
280#endif // __VTK_WRAP__
281
282#define vtkMakeMappedUnstructuredGrid(_className, _impl) \
283 vtkMakeExportedMappedUnstructuredGrid(_className, _impl, )
284
285#define vtkMakeMappedUnstructuredGridWithIter(_className, _impl, _cIter) \
286 vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, )
287
288#endif // vtkMappedUnstructuredGrid_h
289
290// VTK-HeaderTest-Exclude: vtkMappedUnstructuredGrid.h
object to represent cell connectivity
Efficient cell iterator for vtkDataSet topologies.
abstract class to specify cell behavior
Definition vtkCell.h:130
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
virtual vtkCell * GetCell(vtkIdType cellId)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)=0
Topological inquiry to get points defining cell.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:108
Default cell iterator for vtkMappedUnstructuredGrid.
Allows datasets with arbitrary storage layouts to be used with VTK.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
ImplementationType * GetImplementation()
void SetImplementation(ImplementationType *impl)
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Topological inquiry to get points defining cell.
vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[]) override
vtkIdType InternalInsertNextCell(int type, vtkIdList *ptIds) override
int GetMaxCellSize() override
Convenience method returns largest cell size in dataset.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Topological inquiry to get cells using point.
vtkCellIterator * NewCellIterator() override
Return an iterator that traverses the cells in this data set.
vtkSmartPointer< ImplementationType > Impl
void CopyStructure(vtkDataSet *pd) override
Copy the geometric and topological structure of an object.
void InternalReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) override
int IsHomogeneous() override
Traverse cells and determine if cells are all of the same type.
vtkIdType GetNumberOfCells() override
Determine the number of cells composing the dataset.
void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array) override
Fill vtkIdTypeArray container with list of cell Ids.
vtkMappedUnstructuredGrid< Implementation, CellIterator > ThisType
vtkTemplateTypeMacro(SelfType, vtkUnstructuredGridBase)
vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[], vtkCellArray *faces) override
~vtkMappedUnstructuredGrid() override
void Allocate(vtkIdType numCells, int extSize=1000) override
Allocate memory for the number of cells indicated.
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
vtkCell * GetCell(vtkIdType cellId) override
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetCellType(vtkIdType cellId) override
Get type of cell with cellId such that: 0 <= cellId < NumberOfCells.
vtkMTimeType GetMTime() override
Datasets are composite objects and need to check each part for MTime THIS METHOD IS THREAD SAFE.
Allocate and hold a VTK object.
Definition vtkNew.h:167
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types.
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270