VTK  9.0.3
vtkCellIterator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCellIterator.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
15 
65 #ifndef vtkCellIterator_h
66 #define vtkCellIterator_h
67 
68 #include "vtkCellType.h" // For VTK_EMPTY_CELL
69 #include "vtkCommonDataModelModule.h" // For export macro
70 #include "vtkIdList.h" // For inline methods
71 #include "vtkNew.h" // For vtkNew
72 #include "vtkObject.h"
73 
74 class vtkGenericCell;
75 class vtkPoints;
76 
77 class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
78 {
79 public:
80  void PrintSelf(ostream& os, vtkIndent indent) override;
82 
86  void InitTraversal();
87 
91  void GoToNextCell();
92 
96  virtual bool IsDoneWithTraversal() = 0;
97 
102  int GetCellType();
103 
109 
113  virtual vtkIdType GetCellId() = 0;
114 
119  vtkIdList* GetPointIds();
120 
126  vtkPoints* GetPoints();
127 
132  vtkIdList* GetFaces();
133 
139  void GetCell(vtkGenericCell* cell);
140 
145  vtkIdType GetNumberOfPoints();
146 
151  vtkIdType GetNumberOfFaces();
152 
153 protected:
155  ~vtkCellIterator() override;
156 
160  virtual void ResetToFirstCell() = 0;
161 
165  virtual void IncrementToNextCell() = 0;
166 
170  virtual void FetchCellType() = 0;
171 
175  virtual void FetchPointIds() = 0;
176 
180  virtual void FetchPoints() = 0;
181 
188  virtual void FetchFaces() {}
189 
190  int CellType;
194 
195 private:
196  vtkCellIterator(const vtkCellIterator&) = delete;
197  void operator=(const vtkCellIterator&) = delete;
198 
199  enum
200  {
201  UninitializedFlag = 0x0,
202  CellTypeFlag = 0x1,
203  PointIdsFlag = 0x2,
204  PointsFlag = 0x4,
205  FacesFlag = 0x8
206  };
207 
208  void ResetCache()
209  {
210  this->CacheFlags = UninitializedFlag;
211  this->CellType = VTK_EMPTY_CELL;
212  }
213 
214  void SetCache(unsigned char flags) { this->CacheFlags |= flags; }
215 
216  bool CheckCache(unsigned char flags) { return (this->CacheFlags & flags) == flags; }
217 
218  vtkNew<vtkPoints> PointsContainer;
219  vtkNew<vtkIdList> PointIdsContainer;
220  vtkNew<vtkIdList> FacesContainer;
221  unsigned char CacheFlags;
222 };
223 
224 //------------------------------------------------------------------------------
226 {
227  this->ResetToFirstCell();
228  this->ResetCache();
229 }
230 
231 //------------------------------------------------------------------------------
233 {
234  this->IncrementToNextCell();
235  this->ResetCache();
236 }
237 
238 //------------------------------------------------------------------------------
240 {
241  if (!this->CheckCache(CellTypeFlag))
242  {
243  this->FetchCellType();
244  this->SetCache(CellTypeFlag);
245  }
246  return this->CellType;
247 }
248 
249 //------------------------------------------------------------------------------
251 {
252  if (!this->CheckCache(PointIdsFlag))
253  {
254  this->FetchPointIds();
255  this->SetCache(PointIdsFlag);
256  }
257  return this->PointIds;
258 }
259 
260 //------------------------------------------------------------------------------
262 {
263  if (!this->CheckCache(PointsFlag))
264  {
265  this->FetchPoints();
266  this->SetCache(PointsFlag);
267  }
268  return this->Points;
269 }
270 
271 //------------------------------------------------------------------------------
273 {
274  if (!this->CheckCache(FacesFlag))
275  {
276  this->FetchFaces();
277  this->SetCache(FacesFlag);
278  }
279  return this->Faces;
280 }
281 
282 //------------------------------------------------------------------------------
284 {
285  if (!this->CheckCache(PointIdsFlag))
286  {
287  this->FetchPointIds();
288  this->SetCache(PointIdsFlag);
289  }
290  return this->PointIds->GetNumberOfIds();
291 }
292 
293 //------------------------------------------------------------------------------
295 {
296  switch (this->GetCellType())
297  {
298  case VTK_EMPTY_CELL:
299  case VTK_VERTEX:
300  case VTK_POLY_VERTEX:
301  case VTK_LINE:
302  case VTK_POLY_LINE:
303  case VTK_TRIANGLE:
304  case VTK_TRIANGLE_STRIP:
305  case VTK_POLYGON:
306  case VTK_PIXEL:
307  case VTK_QUAD:
308  case VTK_QUADRATIC_EDGE:
310  case VTK_QUADRATIC_QUAD:
315  case VTK_CUBIC_LINE:
325  case VTK_LAGRANGE_CURVE:
328  case VTK_BEZIER_CURVE:
329  case VTK_BEZIER_TRIANGLE:
331  return 0;
332 
333  case VTK_TETRA:
334  case VTK_QUADRATIC_TETRA:
339  return 4;
340 
341  case VTK_PYRAMID:
344  case VTK_WEDGE:
345  case VTK_QUADRATIC_WEDGE:
349  case VTK_LAGRANGE_WEDGE:
350  case VTK_BEZIER_WEDGE:
351  return 5;
352 
353  case VTK_VOXEL:
354  case VTK_HEXAHEDRON:
362  return 6;
363 
365  return 7;
366 
367  case VTK_HEXAGONAL_PRISM:
368  return 8;
369 
370  case VTK_POLYHEDRON: // Need to look these up
371  if (!this->CheckCache(FacesFlag))
372  {
373  this->FetchFaces();
374  this->SetCache(FacesFlag);
375  }
376  return this->Faces->GetNumberOfIds() != 0 ? this->Faces->GetId(0) : 0;
377 
378  default:
379  vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
380  break;
381  }
382 
383  return 0;
384 }
385 
386 #endif // vtkCellIterator_h
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.
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.
vtkIdList * GetFaces()
Get the faces for a polyhedral cell.
void InitTraversal()
Reset to the first cell.
vtkAbstractTypeMacro(vtkCellIterator, vtkObject)
vtkIdList * Faces
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
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.
vtkPoints * Points
~vtkCellIterator() override
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:31
vtkIdType GetNumberOfIds()
Return the number of id's in the list.
Definition: vtkIdList.h:57
vtkIdType GetId(const vtkIdType i)
Return the id at location i.
Definition: vtkIdList.h:62
a simple class to control print indentation
Definition: vtkIndent.h:34
abstract base class for most VTK objects
Definition: vtkObject.h:54
represent and manipulate 3D points
Definition: vtkPoints.h:34
@ VTK_VOXEL
Definition: vtkCellType.h:57
@ VTK_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:70
@ VTK_PARAMETRIC_SURFACE
Definition: vtkCellType.h:92
@ VTK_HIGHER_ORDER_TETRAHEDRON
Definition: vtkCellType.h:103
@ VTK_TRIANGLE_STRIP
Definition: vtkCellType.h:52
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:78
@ VTK_LAGRANGE_CURVE
Definition: vtkCellType.h:109
@ VTK_HIGHER_ORDER_QUAD
Definition: vtkCellType.h:101
@ VTK_PYRAMID
Definition: vtkCellType.h:60
@ VTK_PIXEL
Definition: vtkCellType.h:54
@ VTK_QUADRATIC_WEDGE
Definition: vtkCellType.h:71
@ VTK_BEZIER_WEDGE
Definition: vtkCellType.h:123
@ VTK_BIQUADRATIC_QUAD
Definition: vtkCellType.h:73
@ VTK_HIGHER_ORDER_WEDGE
Definition: vtkCellType.h:104
@ VTK_LAGRANGE_QUADRILATERAL
Definition: vtkCellType.h:111
@ VTK_POLY_LINE
Definition: vtkCellType.h:50
@ VTK_TRIANGLE
Definition: vtkCellType.h:51
@ VTK_BEZIER_TRIANGLE
Definition: vtkCellType.h:119
@ VTK_POLYGON
Definition: vtkCellType.h:53
@ VTK_EMPTY_CELL
Definition: vtkCellType.h:46
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:72
@ VTK_POLYHEDRON
Definition: vtkCellType.h:88
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:74
@ VTK_TETRA
Definition: vtkCellType.h:56
@ VTK_LINE
Definition: vtkCellType.h:49
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:85
@ VTK_BEZIER_HEXAHEDRON
Definition: vtkCellType.h:122
@ VTK_PARAMETRIC_TRI_SURFACE
Definition: vtkCellType.h:93
@ VTK_LAGRANGE_WEDGE
Definition: vtkCellType.h:114
@ VTK_LAGRANGE_HEXAHEDRON
Definition: vtkCellType.h:113
@ VTK_PENTAGONAL_PRISM
Definition: vtkCellType.h:61
@ VTK_HIGHER_ORDER_TRIANGLE
Definition: vtkCellType.h:100
@ VTK_QUADRATIC_QUAD
Definition: vtkCellType.h:67
@ VTK_WEDGE
Definition: vtkCellType.h:59
@ VTK_PARAMETRIC_QUAD_SURFACE
Definition: vtkCellType.h:94
@ VTK_LAGRANGE_TETRAHEDRON
Definition: vtkCellType.h:112
@ VTK_PARAMETRIC_CURVE
Definition: vtkCellType.h:91
@ VTK_BEZIER_CURVE
Definition: vtkCellType.h:118
@ VTK_HIGHER_ORDER_PYRAMID
Definition: vtkCellType.h:105
@ VTK_HEXAGONAL_PRISM
Definition: vtkCellType.h:62
@ VTK_PARAMETRIC_HEX_REGION
Definition: vtkCellType.h:96
@ VTK_BEZIER_QUADRILATERAL
Definition: vtkCellType.h:120
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:76
@ VTK_HEXAHEDRON
Definition: vtkCellType.h:58
@ VTK_CUBIC_LINE
Definition: vtkCellType.h:82
@ VTK_LAGRANGE_TRIANGLE
Definition: vtkCellType.h:110
@ VTK_HIGHER_ORDER_HEXAHEDRON
Definition: vtkCellType.h:106
@ VTK_QUADRATIC_POLYGON
Definition: vtkCellType.h:68
@ VTK_QUAD
Definition: vtkCellType.h:55
@ VTK_QUADRATIC_TRIANGLE
Definition: vtkCellType.h:66
@ VTK_PARAMETRIC_TETRA_REGION
Definition: vtkCellType.h:95
@ VTK_QUADRATIC_EDGE
Definition: vtkCellType.h:65
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:69
@ VTK_HIGHER_ORDER_EDGE
Definition: vtkCellType.h:99
@ VTK_BEZIER_TETRAHEDRON
Definition: vtkCellType.h:121
@ VTK_VERTEX
Definition: vtkCellType.h:47
@ VTK_POLY_VERTEX
Definition: vtkCellType.h:48
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:75
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition: vtkCellType.h:77
@ VTK_HIGHER_ORDER_POLYGON
Definition: vtkCellType.h:102
@ VTK_BIQUADRATIC_TRIANGLE
Definition: vtkCellType.h:79
int vtkIdType
Definition: vtkType.h:338