VTK  9.5.20250619
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 "vtkIdList.h" // For inline methods
162#include "vtkIdTypeArray.h" // For inline methods
163#include "vtkNew.h" // For vtkNew
164#include "vtkObject.h"
165
166VTK_ABI_NAMESPACE_BEGIN
167class vtkGenericCell;
168class vtkPoints;
169
170class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
171{
172public:
173 void PrintSelf(ostream& os, vtkIndent indent) override;
175
179 void InitTraversal();
180
184 void GoToNextCell();
185
189 virtual bool IsDoneWithTraversal() = 0;
190
195 int GetCellType();
196
202
206 virtual vtkIdType GetCellId() = 0;
207
212 vtkIdList* GetPointIds();
213
219 vtkPoints* GetPoints();
220
225 vtkCellArray* GetCellFaces();
226
231 vtkIdList* GetSerializedCellFaces();
232
239
244 vtkIdType GetNumberOfPoints();
245
250 vtkIdType GetNumberOfFaces();
251
252protected:
255
259 virtual void ResetToFirstCell() = 0;
260
264 virtual void IncrementToNextCell() = 0;
265
269 virtual void FetchCellType() = 0;
270
274 virtual void FetchPointIds() = 0;
275
279 virtual void FetchPoints() = 0;
280
287 virtual void FetchFaces() {}
288
293
294private:
295 vtkCellIterator(const vtkCellIterator&) = delete;
296 void operator=(const vtkCellIterator&) = delete;
297
298 enum
299 {
300 UninitializedFlag = 0x0,
301 CellTypeFlag = 0x1,
302 PointIdsFlag = 0x2,
303 PointsFlag = 0x4,
304 FacesFlag = 0x8
305 };
306
307 void ResetCache()
308 {
309 this->CacheFlags = UninitializedFlag;
310 this->CellType = VTK_EMPTY_CELL;
311 }
312
313 void SetCache(unsigned char flags) { this->CacheFlags |= flags; }
314
315 bool CheckCache(unsigned char flags) { return (this->CacheFlags & flags) == flags; }
316
317 vtkNew<vtkPoints> PointsContainer;
318 vtkNew<vtkIdList> PointIdsContainer;
319 vtkNew<vtkCellArray> FacesContainer;
320 vtkNew<vtkIdList> LegacyFacesContainer;
321 unsigned char CacheFlags;
322};
323
324//------------------------------------------------------------------------------
326{
327 this->ResetToFirstCell();
328 this->ResetCache();
329}
330
331//------------------------------------------------------------------------------
333{
334 this->IncrementToNextCell();
335 this->ResetCache();
336}
337
338//------------------------------------------------------------------------------
340{
341 if (!this->CheckCache(CellTypeFlag))
342 {
343 this->FetchCellType();
344 this->SetCache(CellTypeFlag);
345 }
346 return this->CellType;
347}
348
349//------------------------------------------------------------------------------
351{
352 if (!this->CheckCache(PointIdsFlag))
353 {
354 this->FetchPointIds();
355 this->SetCache(PointIdsFlag);
356 }
357 return this->PointIds;
358}
359
360//------------------------------------------------------------------------------
362{
363 if (!this->CheckCache(PointsFlag))
364 {
365 this->FetchPoints();
366 this->SetCache(PointsFlag);
367 }
368 return this->Points;
369}
370
371//------------------------------------------------------------------------------
373{
374 if (!this->CheckCache(FacesFlag))
375 {
376 this->FetchFaces();
377 this->SetCache(FacesFlag);
378 }
379 return this->Faces;
380}
381
382//------------------------------------------------------------------------------
384{
385 if (!this->CheckCache(FacesFlag))
386 {
387 this->FetchFaces();
388 this->SetCache(FacesFlag);
389 }
390 // Export Legacy Format
392 this->Faces->ExportLegacyFormat(tmp);
393 this->LegacyFacesContainer->Initialize();
394 this->LegacyFacesContainer->InsertNextId(this->Faces->GetNumberOfCells());
395 for (vtkIdType idx = 0; idx < tmp->GetNumberOfValues(); ++idx)
396 {
397 this->LegacyFacesContainer->InsertNextId(tmp->GetValue(idx));
398 }
399 return this->LegacyFacesContainer;
400}
401
402//------------------------------------------------------------------------------
404{
405 if (!this->CheckCache(PointIdsFlag))
406 {
407 this->FetchPointIds();
408 this->SetCache(PointIdsFlag);
409 }
410 return this->PointIds->GetNumberOfIds();
411}
412
413//------------------------------------------------------------------------------
415{
416 switch (this->GetCellType())
417 {
418 case VTK_EMPTY_CELL:
419 case VTK_VERTEX:
420 case VTK_POLY_VERTEX:
421 case VTK_LINE:
422 case VTK_POLY_LINE:
423 case VTK_TRIANGLE:
425 case VTK_POLYGON:
426 case VTK_PIXEL:
427 case VTK_QUAD:
435 case VTK_CUBIC_LINE:
448 case VTK_BEZIER_CURVE:
451 return 0;
452
453 case VTK_TETRA:
459 return 4;
460
461 case VTK_PYRAMID:
465 case VTK_WEDGE:
471 case VTK_BEZIER_WEDGE:
472 return 5;
473
474 case VTK_VOXEL:
475 case VTK_HEXAHEDRON:
483 return 6;
484
486 return 7;
487
489 return 8;
490
491 case VTK_POLYHEDRON: // Need to look these up
492 if (!this->CheckCache(FacesFlag))
493 {
494 this->FetchFaces();
495 this->SetCache(FacesFlag);
496 }
497 return this->Faces->GetNumberOfCells();
498
499 default:
500 vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
501 break;
502 }
503
504 return 0;
505}
506
507VTK_ABI_NAMESPACE_END
508#endif // vtkCellIterator_h
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.
void GetCell(vtkGenericCell *cell)
Write the current full cell information into the argument.
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.
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.
~vtkCellIterator() override
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition vtkIdList.h:159
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:336
a simple class to control print indentation
Definition vtkIndent.h:108
Allocate and hold a VTK object.
Definition vtkNew.h:167
abstract base class for most VTK objects
Definition vtkObject.h:162
represent and manipulate 3D points
Definition vtkPoints.h:139
@ VTK_VOXEL
Definition vtkCellType.h:48
@ VTK_QUADRATIC_HEXAHEDRON
Definition vtkCellType.h:61
@ VTK_PARAMETRIC_SURFACE
Definition vtkCellType.h:84
@ VTK_HIGHER_ORDER_TETRAHEDRON
Definition vtkCellType.h:95
@ VTK_TRIANGLE_STRIP
Definition vtkCellType.h:43
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition vtkCellType.h:70
@ VTK_LAGRANGE_CURVE
@ VTK_HIGHER_ORDER_QUAD
Definition vtkCellType.h:93
@ VTK_PYRAMID
Definition vtkCellType.h:51
@ VTK_PIXEL
Definition vtkCellType.h:45
@ VTK_QUADRATIC_WEDGE
Definition vtkCellType.h:62
@ VTK_BEZIER_WEDGE
@ VTK_BIQUADRATIC_QUAD
Definition vtkCellType.h:64
@ VTK_HIGHER_ORDER_WEDGE
Definition vtkCellType.h:96
@ VTK_LAGRANGE_QUADRILATERAL
@ VTK_POLY_LINE
Definition vtkCellType.h:41
@ VTK_TRIQUADRATIC_PYRAMID
Definition vtkCellType.h:66
@ VTK_TRIANGLE
Definition vtkCellType.h:42
@ VTK_BEZIER_TRIANGLE
@ VTK_POLYGON
Definition vtkCellType.h:44
@ VTK_EMPTY_CELL
Definition vtkCellType.h:37
@ VTK_QUADRATIC_PYRAMID
Definition vtkCellType.h:63
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition vtkCellType.h:65
@ VTK_TETRA
Definition vtkCellType.h:47
@ VTK_LINE
Definition vtkCellType.h:40
@ VTK_CONVEX_POINT_SET
Definition vtkCellType.h:77
@ VTK_BEZIER_HEXAHEDRON
@ VTK_PARAMETRIC_TRI_SURFACE
Definition vtkCellType.h:85
@ VTK_LAGRANGE_WEDGE
@ VTK_LAGRANGE_HEXAHEDRON
@ VTK_PENTAGONAL_PRISM
Definition vtkCellType.h:52
@ VTK_HIGHER_ORDER_TRIANGLE
Definition vtkCellType.h:92
@ VTK_QUADRATIC_QUAD
Definition vtkCellType.h:58
@ VTK_WEDGE
Definition vtkCellType.h:50
@ VTK_PARAMETRIC_QUAD_SURFACE
Definition vtkCellType.h:86
@ VTK_LAGRANGE_TETRAHEDRON
@ VTK_PARAMETRIC_CURVE
Definition vtkCellType.h:83
@ VTK_BEZIER_CURVE
@ VTK_HIGHER_ORDER_PYRAMID
Definition vtkCellType.h:97
@ VTK_HEXAGONAL_PRISM
Definition vtkCellType.h:53
@ VTK_PARAMETRIC_HEX_REGION
Definition vtkCellType.h:88
@ VTK_BEZIER_QUADRILATERAL
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition vtkCellType.h:68
@ VTK_HEXAHEDRON
Definition vtkCellType.h:49
@ VTK_CUBIC_LINE
Definition vtkCellType.h:74
@ VTK_LAGRANGE_TRIANGLE
@ VTK_HIGHER_ORDER_HEXAHEDRON
Definition vtkCellType.h:98
@ VTK_QUADRATIC_POLYGON
Definition vtkCellType.h:59
@ VTK_QUAD
Definition vtkCellType.h:46
@ VTK_QUADRATIC_TRIANGLE
Definition vtkCellType.h:57
@ VTK_PARAMETRIC_TETRA_REGION
Definition vtkCellType.h:87
@ VTK_QUADRATIC_EDGE
Definition vtkCellType.h:56
@ VTK_QUADRATIC_TETRA
Definition vtkCellType.h:60
@ VTK_HIGHER_ORDER_EDGE
Definition vtkCellType.h:91
@ VTK_BEZIER_TETRAHEDRON
@ VTK_VERTEX
Definition vtkCellType.h:38
@ VTK_POLY_VERTEX
Definition vtkCellType.h:39
@ VTK_QUADRATIC_LINEAR_QUAD
Definition vtkCellType.h:67
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition vtkCellType.h:69
@ VTK_HIGHER_ORDER_POLYGON
Definition vtkCellType.h:94
@ VTK_BIQUADRATIC_TRIANGLE
Definition vtkCellType.h:71
int vtkIdType
Definition vtkType.h:332