VTK  9.3.20240327
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 
155 #ifndef vtkCellIterator_h
156 #define vtkCellIterator_h
157 
158 #include "vtkCellArray.h" // For inline methods
159 #include "vtkCellType.h" // For VTK_EMPTY_CELL
160 #include "vtkCommonDataModelModule.h" // For export macro
161 #include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_4_0
162 #include "vtkIdList.h" // For inline methods
163 #include "vtkIdTypeArray.h" // For inline methods
164 #include "vtkNew.h" // For vtkNew
165 #include "vtkObject.h"
166 
167 VTK_ABI_NAMESPACE_BEGIN
168 class vtkGenericCell;
169 class vtkPoints;
170 
171 class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
172 {
173 public:
174  void PrintSelf(ostream& os, vtkIndent indent) override;
176 
180  void InitTraversal();
181 
185  void GoToNextCell();
186 
190  virtual bool IsDoneWithTraversal() = 0;
191 
196  int GetCellType();
197 
203 
207  virtual vtkIdType GetCellId() = 0;
208 
213  vtkIdList* GetPointIds();
214 
220  vtkPoints* GetPoints();
221 
226  vtkCellArray* GetCellFaces();
227 
232  VTK_DEPRECATED_IN_9_4_0("Please use GetCellFaces instead.")
233  vtkIdList* GetFaces();
234 
240  void GetCell(vtkGenericCell* cell);
241 
246  vtkIdType GetNumberOfPoints();
247 
252  vtkIdType GetNumberOfFaces();
253 
254 protected:
256  ~vtkCellIterator() override;
257 
261  virtual void ResetToFirstCell() = 0;
262 
266  virtual void IncrementToNextCell() = 0;
267 
271  virtual void FetchCellType() = 0;
272 
276  virtual void FetchPointIds() = 0;
277 
281  virtual void FetchPoints() = 0;
282 
289  virtual void FetchFaces() {}
290 
291  int CellType;
295 
296 private:
297  vtkCellIterator(const vtkCellIterator&) = delete;
298  void operator=(const vtkCellIterator&) = delete;
299 
300  enum
301  {
302  UninitializedFlag = 0x0,
303  CellTypeFlag = 0x1,
304  PointIdsFlag = 0x2,
305  PointsFlag = 0x4,
306  FacesFlag = 0x8
307  };
308 
309  void ResetCache()
310  {
311  this->CacheFlags = UninitializedFlag;
312  this->CellType = VTK_EMPTY_CELL;
313  }
314 
315  void SetCache(unsigned char flags) { this->CacheFlags |= flags; }
316 
317  bool CheckCache(unsigned char flags) { return (this->CacheFlags & flags) == flags; }
318 
319  vtkNew<vtkPoints> PointsContainer;
320  vtkNew<vtkIdList> PointIdsContainer;
321  vtkNew<vtkCellArray> FacesContainer;
322  vtkNew<vtkIdList> LegacyFacesContainer;
323  unsigned char CacheFlags;
324 };
325 
326 //------------------------------------------------------------------------------
328 {
329  this->ResetToFirstCell();
330  this->ResetCache();
331 }
332 
333 //------------------------------------------------------------------------------
335 {
336  this->IncrementToNextCell();
337  this->ResetCache();
338 }
339 
340 //------------------------------------------------------------------------------
342 {
343  if (!this->CheckCache(CellTypeFlag))
344  {
345  this->FetchCellType();
346  this->SetCache(CellTypeFlag);
347  }
348  return this->CellType;
349 }
350 
351 //------------------------------------------------------------------------------
353 {
354  if (!this->CheckCache(PointIdsFlag))
355  {
356  this->FetchPointIds();
357  this->SetCache(PointIdsFlag);
358  }
359  return this->PointIds;
360 }
361 
362 //------------------------------------------------------------------------------
364 {
365  if (!this->CheckCache(PointsFlag))
366  {
367  this->FetchPoints();
368  this->SetCache(PointsFlag);
369  }
370  return this->Points;
371 }
372 
373 //------------------------------------------------------------------------------
375 {
376  if (!this->CheckCache(FacesFlag))
377  {
378  this->FetchFaces();
379  this->SetCache(FacesFlag);
380  }
381  return this->Faces;
382 }
383 
384 //------------------------------------------------------------------------------
385 // To be removed when deprecating
387 {
388  if (!this->CheckCache(FacesFlag))
389  {
390  this->FetchFaces();
391  this->SetCache(FacesFlag);
392  }
393  // Export Legacy Format
395  this->Faces->ExportLegacyFormat(tmp);
396  this->LegacyFacesContainer->Initialize();
397  this->LegacyFacesContainer->InsertNextId(this->Faces->GetNumberOfCells());
398  for (vtkIdType idx = 0; idx < tmp->GetNumberOfValues(); ++idx)
399  {
400  this->LegacyFacesContainer->InsertNextId(tmp->GetValue(idx));
401  }
402  return this->LegacyFacesContainer;
403 }
404 
405 //------------------------------------------------------------------------------
407 {
408  if (!this->CheckCache(PointIdsFlag))
409  {
410  this->FetchPointIds();
411  this->SetCache(PointIdsFlag);
412  }
413  return this->PointIds->GetNumberOfIds();
414 }
415 
416 //------------------------------------------------------------------------------
418 {
419  switch (this->GetCellType())
420  {
421  case VTK_EMPTY_CELL:
422  case VTK_VERTEX:
423  case VTK_POLY_VERTEX:
424  case VTK_LINE:
425  case VTK_POLY_LINE:
426  case VTK_TRIANGLE:
427  case VTK_TRIANGLE_STRIP:
428  case VTK_POLYGON:
429  case VTK_PIXEL:
430  case VTK_QUAD:
431  case VTK_QUADRATIC_EDGE:
433  case VTK_QUADRATIC_QUAD:
438  case VTK_CUBIC_LINE:
448  case VTK_LAGRANGE_CURVE:
451  case VTK_BEZIER_CURVE:
452  case VTK_BEZIER_TRIANGLE:
454  return 0;
455 
456  case VTK_TETRA:
457  case VTK_QUADRATIC_TETRA:
462  return 4;
463 
464  case VTK_PYRAMID:
468  case VTK_WEDGE:
469  case VTK_QUADRATIC_WEDGE:
473  case VTK_LAGRANGE_WEDGE:
474  case VTK_BEZIER_WEDGE:
475  return 5;
476 
477  case VTK_VOXEL:
478  case VTK_HEXAHEDRON:
486  return 6;
487 
489  return 7;
490 
491  case VTK_HEXAGONAL_PRISM:
492  return 8;
493 
494  case VTK_POLYHEDRON: // Need to look these up
495  if (!this->CheckCache(FacesFlag))
496  {
497  this->FetchFaces();
498  this->SetCache(FacesFlag);
499  }
500  return this->Faces->GetNumberOfCells();
501 
502  default:
503  vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
504  break;
505  }
506 
507  return 0;
508 }
509 
510 VTK_ABI_NAMESPACE_END
511 #endif // vtkCellIterator_h
ValueType GetValue(vtkIdType valueIdx) const
Get the value at valueIdx.
vtkIdType GetNumberOfValues() const
Get the total number of values in the array.
object to represent cell connectivity
Definition: vtkCellArray.h:285
vtkIdType GetNumberOfCells() const override
Get the number of cells in the array.
Definition: vtkCellArray.h:418
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).
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.
vtkIdList * GetFaces()
Get the faces for a polyhedral 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.
vtkIdList * PointIds
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.
vtkCellArray * Faces
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.
vtkPoints * Points
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:132
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition: vtkIdList.h:158
void Initialize()
Release memory and restore to unallocated state.
vtkIdType InsertNextId(vtkIdType vtkid)
Add the id specified to the end of the list.
Definition: vtkIdList.h:335
a simple class to control print indentation
Definition: vtkIndent.h:108
abstract base class for most VTK objects
Definition: vtkObject.h:161
represent and manipulate 3D points
Definition: vtkPoints.h:138
VTKIOCATALYSTCONDUIT_EXPORT int GetCellType(const std::string &shape)
Get vtk cell type from conduit shape name throw a runtime_error on unsupported type.
@ VTK_VOXEL
Definition: vtkCellType.h:67
@ VTK_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:80
@ VTK_PARAMETRIC_SURFACE
Definition: vtkCellType.h:103
@ VTK_HIGHER_ORDER_TETRAHEDRON
Definition: vtkCellType.h:114
@ VTK_TRIANGLE_STRIP
Definition: vtkCellType.h:62
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:89
@ VTK_LAGRANGE_CURVE
Definition: vtkCellType.h:120
@ VTK_HIGHER_ORDER_QUAD
Definition: vtkCellType.h:112
@ VTK_PYRAMID
Definition: vtkCellType.h:70
@ VTK_PIXEL
Definition: vtkCellType.h:64
@ VTK_QUADRATIC_WEDGE
Definition: vtkCellType.h:81
@ VTK_BEZIER_WEDGE
Definition: vtkCellType.h:134
@ VTK_BIQUADRATIC_QUAD
Definition: vtkCellType.h:83
@ VTK_HIGHER_ORDER_WEDGE
Definition: vtkCellType.h:115
@ VTK_LAGRANGE_QUADRILATERAL
Definition: vtkCellType.h:122
@ VTK_POLY_LINE
Definition: vtkCellType.h:60
@ VTK_TRIQUADRATIC_PYRAMID
Definition: vtkCellType.h:85
@ VTK_TRIANGLE
Definition: vtkCellType.h:61
@ VTK_BEZIER_TRIANGLE
Definition: vtkCellType.h:130
@ VTK_POLYGON
Definition: vtkCellType.h:63
@ VTK_EMPTY_CELL
Definition: vtkCellType.h:56
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:82
@ VTK_POLYHEDRON
Definition: vtkCellType.h:99
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:84
@ VTK_TETRA
Definition: vtkCellType.h:66
@ VTK_LINE
Definition: vtkCellType.h:59
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:96
@ VTK_BEZIER_HEXAHEDRON
Definition: vtkCellType.h:133
@ VTK_PARAMETRIC_TRI_SURFACE
Definition: vtkCellType.h:104
@ VTK_LAGRANGE_WEDGE
Definition: vtkCellType.h:125
@ VTK_LAGRANGE_HEXAHEDRON
Definition: vtkCellType.h:124
@ VTK_PENTAGONAL_PRISM
Definition: vtkCellType.h:71
@ VTK_HIGHER_ORDER_TRIANGLE
Definition: vtkCellType.h:111
@ VTK_QUADRATIC_QUAD
Definition: vtkCellType.h:77
@ VTK_WEDGE
Definition: vtkCellType.h:69
@ VTK_PARAMETRIC_QUAD_SURFACE
Definition: vtkCellType.h:105
@ VTK_LAGRANGE_TETRAHEDRON
Definition: vtkCellType.h:123
@ VTK_PARAMETRIC_CURVE
Definition: vtkCellType.h:102
@ VTK_BEZIER_CURVE
Definition: vtkCellType.h:129
@ VTK_HIGHER_ORDER_PYRAMID
Definition: vtkCellType.h:116
@ VTK_HEXAGONAL_PRISM
Definition: vtkCellType.h:72
@ VTK_PARAMETRIC_HEX_REGION
Definition: vtkCellType.h:107
@ VTK_BEZIER_QUADRILATERAL
Definition: vtkCellType.h:131
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:87
@ VTK_HEXAHEDRON
Definition: vtkCellType.h:68
@ VTK_CUBIC_LINE
Definition: vtkCellType.h:93
@ VTK_LAGRANGE_TRIANGLE
Definition: vtkCellType.h:121
@ VTK_HIGHER_ORDER_HEXAHEDRON
Definition: vtkCellType.h:117
@ VTK_QUADRATIC_POLYGON
Definition: vtkCellType.h:78
@ VTK_QUAD
Definition: vtkCellType.h:65
@ VTK_QUADRATIC_TRIANGLE
Definition: vtkCellType.h:76
@ VTK_PARAMETRIC_TETRA_REGION
Definition: vtkCellType.h:106
@ VTK_QUADRATIC_EDGE
Definition: vtkCellType.h:75
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:79
@ VTK_HIGHER_ORDER_EDGE
Definition: vtkCellType.h:110
@ VTK_BEZIER_TETRAHEDRON
Definition: vtkCellType.h:132
@ VTK_VERTEX
Definition: vtkCellType.h:57
@ VTK_POLY_VERTEX
Definition: vtkCellType.h:58
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:86
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition: vtkCellType.h:88
@ VTK_HIGHER_ORDER_POLYGON
Definition: vtkCellType.h:113
@ VTK_BIQUADRATIC_TRIANGLE
Definition: vtkCellType.h:90
#define VTK_DEPRECATED_IN_9_4_0(reason)
int vtkIdType
Definition: vtkType.h:315