VTK  9.3.1
vtkCellIterator.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
3 
56 #ifndef vtkCellIterator_h
57 #define vtkCellIterator_h
58 
59 #include "vtkCellType.h" // For VTK_EMPTY_CELL
60 #include "vtkCommonDataModelModule.h" // For export macro
61 #include "vtkIdList.h" // For inline methods
62 #include "vtkNew.h" // For vtkNew
63 #include "vtkObject.h"
64 
65 VTK_ABI_NAMESPACE_BEGIN
66 class vtkGenericCell;
67 class vtkPoints;
68 
69 class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
70 {
71 public:
72  void PrintSelf(ostream& os, vtkIndent indent) override;
73  vtkAbstractTypeMacro(vtkCellIterator, vtkObject);
74 
78  void InitTraversal();
79 
83  void GoToNextCell();
84 
88  virtual bool IsDoneWithTraversal() = 0;
89 
94  int GetCellType();
95 
100  int GetCellDimension();
101 
105  virtual vtkIdType GetCellId() = 0;
106 
111  vtkIdList* GetPointIds();
112 
118  vtkPoints* GetPoints();
119 
124  vtkIdList* GetFaces();
125 
131  void GetCell(vtkGenericCell* cell);
132 
137  vtkIdType GetNumberOfPoints();
138 
143  vtkIdType GetNumberOfFaces();
144 
145 protected:
146  vtkCellIterator();
147  ~vtkCellIterator() override;
148 
152  virtual void ResetToFirstCell() = 0;
153 
157  virtual void IncrementToNextCell() = 0;
158 
162  virtual void FetchCellType() = 0;
163 
167  virtual void FetchPointIds() = 0;
168 
172  virtual void FetchPoints() = 0;
173 
180  virtual void FetchFaces() {}
181 
182  int CellType;
186 
187 private:
188  vtkCellIterator(const vtkCellIterator&) = delete;
189  void operator=(const vtkCellIterator&) = delete;
190 
191  enum
192  {
193  UninitializedFlag = 0x0,
194  CellTypeFlag = 0x1,
195  PointIdsFlag = 0x2,
196  PointsFlag = 0x4,
197  FacesFlag = 0x8
198  };
199 
200  void ResetCache()
201  {
202  this->CacheFlags = UninitializedFlag;
203  this->CellType = VTK_EMPTY_CELL;
204  }
205 
206  void SetCache(unsigned char flags) { this->CacheFlags |= flags; }
207 
208  bool CheckCache(unsigned char flags) { return (this->CacheFlags & flags) == flags; }
209 
210  vtkNew<vtkPoints> PointsContainer;
211  vtkNew<vtkIdList> PointIdsContainer;
212  vtkNew<vtkIdList> FacesContainer;
213  unsigned char CacheFlags;
214 };
215 
216 //------------------------------------------------------------------------------
218 {
219  this->ResetToFirstCell();
220  this->ResetCache();
221 }
222 
223 //------------------------------------------------------------------------------
225 {
226  this->IncrementToNextCell();
227  this->ResetCache();
228 }
229 
230 //------------------------------------------------------------------------------
232 {
233  if (!this->CheckCache(CellTypeFlag))
234  {
235  this->FetchCellType();
236  this->SetCache(CellTypeFlag);
237  }
238  return this->CellType;
239 }
240 
241 //------------------------------------------------------------------------------
243 {
244  if (!this->CheckCache(PointIdsFlag))
245  {
246  this->FetchPointIds();
247  this->SetCache(PointIdsFlag);
248  }
249  return this->PointIds;
250 }
251 
252 //------------------------------------------------------------------------------
254 {
255  if (!this->CheckCache(PointsFlag))
256  {
257  this->FetchPoints();
258  this->SetCache(PointsFlag);
259  }
260  return this->Points;
261 }
262 
263 //------------------------------------------------------------------------------
265 {
266  if (!this->CheckCache(FacesFlag))
267  {
268  this->FetchFaces();
269  this->SetCache(FacesFlag);
270  }
271  return this->Faces;
272 }
273 
274 //------------------------------------------------------------------------------
276 {
277  if (!this->CheckCache(PointIdsFlag))
278  {
279  this->FetchPointIds();
280  this->SetCache(PointIdsFlag);
281  }
282  return this->PointIds->GetNumberOfIds();
283 }
284 
285 //------------------------------------------------------------------------------
287 {
288  switch (this->GetCellType())
289  {
290  case VTK_EMPTY_CELL:
291  case VTK_VERTEX:
292  case VTK_POLY_VERTEX:
293  case VTK_LINE:
294  case VTK_POLY_LINE:
295  case VTK_TRIANGLE:
296  case VTK_TRIANGLE_STRIP:
297  case VTK_POLYGON:
298  case VTK_PIXEL:
299  case VTK_QUAD:
300  case VTK_QUADRATIC_EDGE:
302  case VTK_QUADRATIC_QUAD:
307  case VTK_CUBIC_LINE:
317  case VTK_LAGRANGE_CURVE:
320  case VTK_BEZIER_CURVE:
321  case VTK_BEZIER_TRIANGLE:
323  return 0;
324 
325  case VTK_TETRA:
326  case VTK_QUADRATIC_TETRA:
331  return 4;
332 
333  case VTK_PYRAMID:
337  case VTK_WEDGE:
338  case VTK_QUADRATIC_WEDGE:
342  case VTK_LAGRANGE_WEDGE:
343  case VTK_BEZIER_WEDGE:
344  return 5;
345 
346  case VTK_VOXEL:
347  case VTK_HEXAHEDRON:
355  return 6;
356 
358  return 7;
359 
360  case VTK_HEXAGONAL_PRISM:
361  return 8;
362 
363  case VTK_POLYHEDRON: // Need to look these up
364  if (!this->CheckCache(FacesFlag))
365  {
366  this->FetchFaces();
367  this->SetCache(FacesFlag);
368  }
369  return this->Faces->GetNumberOfIds() != 0 ? this->Faces->GetId(0) : 0;
370 
371  default:
372  vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
373  break;
374  }
375 
376  return 0;
377 }
378 
379 VTK_ABI_NAMESPACE_END
380 #endif // vtkCellIterator_h
vtkIdType GetNumberOfFaces()
Return the number of faces in the current cell.
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.
vtkPoints * Points
int vtkIdType
Definition: vtkType.h:315
virtual void FetchCellType()=0
Lookup the cell type in the data set and store it in this->CellType.
virtual void ResetToFirstCell()=0
Update internal state to point to the first cell.
provides thread-safe access to cells
virtual void IncrementToNextCell()=0
Update internal state to point to the next cell.
vtkIdType GetId(vtkIdType i)
Return the id at location i.
Definition: vtkIdList.h:54
void GoToNextCell()
Increment to next cell.
virtual void FetchPointIds()=0
Lookup the cell point ids in the data set and store them in this->PointIds.
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
vtkIdList * GetPointIds()
Get the ids of the points in the current cell.
vtkIdList * GetFaces()
Get the faces for a polyhedral cell.
virtual void FetchFaces()
Lookup the cell faces in the data set and store them in this->Faces.
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition: vtkIdList.h:49
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
int GetCellType()
Get the current cell type (e.g.
void InitTraversal()
Reset to the first cell.
vtkPoints * GetPoints()
Get the points in the current cell.
virtual void FetchPoints()=0
Lookup the cell points in the data set and store them in this->Points.
vtkIdList * Faces
Efficient cell iterator for vtkDataSet topologies.
vtkIdList * PointIds
represent and manipulate 3D points
Definition: vtkPoints.h:28
vtkIdType GetNumberOfPoints()
Return the number of points in the current cell.