155#ifndef vtkCellIterator_h
156#define vtkCellIterator_h
160#include "vtkCommonDataModelModule.h"
167VTK_ABI_NAMESPACE_BEGIN
180 void InitTraversal();
267 virtual
void ResetToFirstCell() = 0;
272 virtual
void IncrementToNextCell() = 0;
277 virtual
void FetchCellType() = 0;
282 virtual
void FetchPointIds() = 0;
287 virtual
void FetchPoints() = 0;
295 virtual
void FetchFaces() {}
308 UninitializedFlag = 0x0,
317 this->CacheFlags = UninitializedFlag;
321 void SetCache(
unsigned char flags) { this->CacheFlags |= flags; }
323 bool CheckCache(
unsigned char flags) {
return (this->CacheFlags & flags) == flags; }
329 unsigned char CacheFlags;
349 if (!this->CheckCache(CellTypeFlag))
352 this->SetCache(CellTypeFlag);
360 if (!this->CheckCache(PointIdsFlag))
363 this->SetCache(PointIdsFlag);
371 if (!this->CheckCache(PointsFlag))
374 this->SetCache(PointsFlag);
382 if (!this->CheckCache(FacesFlag))
385 this->SetCache(FacesFlag);
393 if (!this->CheckCache(FacesFlag))
396 this->SetCache(FacesFlag);
403 for (
vtkIdType idx = 0; idx < tmp->GetNumberOfValues(); ++idx)
405 this->LegacyFacesContainer->
InsertNextId(tmp->GetValue(idx));
407 return this->LegacyFacesContainer;
413 if (!this->CheckCache(PointIdsFlag))
416 this->SetCache(PointIdsFlag);
500 if (!this->CheckCache(FacesFlag))
503 this->SetCache(FacesFlag);
508 vtkGenericWarningMacro(
"Unknown cell type: " << this->CellType);
object to represent cell connectivity
vtkIdType GetNumberOfCells() const override
Get the number of cells in the array.
void ExportLegacyFormat(vtkIdTypeArray *data)
Fill data with the old-style vtkCellArray data layout, e.g.
Efficient cell iterator for vtkDataSet topologies.
int GetCellDimension()
Get the current cell dimension (0, 1, 2, or 3).
vtkIdList * GetSerializedCellFaces()
Get a serialized view of the faces for a polyhedral cell.
virtual void FetchPoints()=0
Lookup the cell points in the data set and store them in this->Points.
virtual void FetchPointIds()=0
Lookup the cell point ids in the data set and store them in this->PointIds.
vtkIdType GetNumberOfFaces()
Return the number of faces in the current cell.
void InitTraversal()
Reset to the first cell.
vtkAbstractTypeMacro(vtkCellIterator, vtkObject)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkPoints * GetPoints()
Get the points in the current cell.
vtkCellArray * GetCellFaces()
Get the faces for a polyhedral cell.
virtual void FetchCellType()=0
Lookup the cell type in the data set and store it in this->CellType.
vtkIdType GetNumberOfPoints()
Return the number of points in the current cell.
virtual void IncrementToNextCell()=0
Update internal state to point to the next cell.
virtual vtkIdType GetCellId()=0
Get the id of the current cell.
virtual void ResetToFirstCell()=0
Update internal state to point to the first cell.
vtkIdList * GetPointIds()
Get the ids of the points in the current cell.
int GetCellType()
Get the current cell type (e.g.
virtual void FetchFaces()
Lookup the cell faces in the data set and store them in this->Faces.
void GoToNextCell()
Increment to next cell.
virtual bool IsDoneWithTraversal()=0
Returns false while the iterator is valid.
provides thread-safe access to cells
list of point or cell ids
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
void Initialize()
Release memory and restore to unallocated state.
vtkIdType InsertNextId(vtkIdType vtkid)
Add the id specified to the end of the list.
a simple class to control print indentation
Allocate and hold a VTK object.
abstract base class for most VTK objects
represent and manipulate 3D points
@ VTK_QUADRATIC_HEXAHEDRON
@ VTK_HIGHER_ORDER_TETRAHEDRON
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
@ VTK_LAGRANGE_QUADRILATERAL
@ VTK_TRIQUADRATIC_PYRAMID
@ VTK_TRIQUADRATIC_HEXAHEDRON
@ VTK_PARAMETRIC_TRI_SURFACE
@ VTK_LAGRANGE_HEXAHEDRON
@ VTK_HIGHER_ORDER_TRIANGLE
@ VTK_PARAMETRIC_QUAD_SURFACE
@ VTK_LAGRANGE_TETRAHEDRON
@ VTK_HIGHER_ORDER_PYRAMID
@ VTK_PARAMETRIC_HEX_REGION
@ VTK_BEZIER_QUADRILATERAL
@ VTK_QUADRATIC_LINEAR_WEDGE
@ VTK_HIGHER_ORDER_HEXAHEDRON
@ VTK_PARAMETRIC_TETRA_REGION
@ VTK_QUADRATIC_LINEAR_QUAD
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
@ VTK_HIGHER_ORDER_POLYGON
@ VTK_BIQUADRATIC_TRIANGLE
#define VTK_DEPRECATED_IN_9_4_0(reason)