VTK  9.1.0
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 
164 #ifndef vtkCellIterator_h
165 #define vtkCellIterator_h
166 
167 #include "vtkCellType.h" // For VTK_EMPTY_CELL
168 #include "vtkCommonDataModelModule.h" // For export macro
169 #include "vtkIdList.h" // For inline methods
170 #include "vtkNew.h" // For vtkNew
171 #include "vtkObject.h"
172 
173 class vtkGenericCell;
174 class vtkPoints;
175 
176 class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
177 {
178 public:
179  void PrintSelf(ostream& os, vtkIndent indent) override;
181 
185  void InitTraversal();
186 
190  void GoToNextCell();
191 
195  virtual bool IsDoneWithTraversal() = 0;
196 
201  int GetCellType();
202 
208 
212  virtual vtkIdType GetCellId() = 0;
213 
218  vtkIdList* GetPointIds();
219 
225  vtkPoints* GetPoints();
226 
231  vtkIdList* GetFaces();
232 
238  void GetCell(vtkGenericCell* cell);
239 
244  vtkIdType GetNumberOfPoints();
245 
250  vtkIdType GetNumberOfFaces();
251 
252 protected:
254  ~vtkCellIterator() override;
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 
289  int CellType;
293 
294 private:
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<vtkIdList> FacesContainer;
320  unsigned char CacheFlags;
321 };
322 
323 //------------------------------------------------------------------------------
325 {
326  this->ResetToFirstCell();
327  this->ResetCache();
328 }
329 
330 //------------------------------------------------------------------------------
332 {
333  this->IncrementToNextCell();
334  this->ResetCache();
335 }
336 
337 //------------------------------------------------------------------------------
339 {
340  if (!this->CheckCache(CellTypeFlag))
341  {
342  this->FetchCellType();
343  this->SetCache(CellTypeFlag);
344  }
345  return this->CellType;
346 }
347 
348 //------------------------------------------------------------------------------
350 {
351  if (!this->CheckCache(PointIdsFlag))
352  {
353  this->FetchPointIds();
354  this->SetCache(PointIdsFlag);
355  }
356  return this->PointIds;
357 }
358 
359 //------------------------------------------------------------------------------
361 {
362  if (!this->CheckCache(PointsFlag))
363  {
364  this->FetchPoints();
365  this->SetCache(PointsFlag);
366  }
367  return this->Points;
368 }
369 
370 //------------------------------------------------------------------------------
372 {
373  if (!this->CheckCache(FacesFlag))
374  {
375  this->FetchFaces();
376  this->SetCache(FacesFlag);
377  }
378  return this->Faces;
379 }
380 
381 //------------------------------------------------------------------------------
383 {
384  if (!this->CheckCache(PointIdsFlag))
385  {
386  this->FetchPointIds();
387  this->SetCache(PointIdsFlag);
388  }
389  return this->PointIds->GetNumberOfIds();
390 }
391 
392 //------------------------------------------------------------------------------
394 {
395  switch (this->GetCellType())
396  {
397  case VTK_EMPTY_CELL:
398  case VTK_VERTEX:
399  case VTK_POLY_VERTEX:
400  case VTK_LINE:
401  case VTK_POLY_LINE:
402  case VTK_TRIANGLE:
403  case VTK_TRIANGLE_STRIP:
404  case VTK_POLYGON:
405  case VTK_PIXEL:
406  case VTK_QUAD:
407  case VTK_QUADRATIC_EDGE:
409  case VTK_QUADRATIC_QUAD:
414  case VTK_CUBIC_LINE:
424  case VTK_LAGRANGE_CURVE:
427  case VTK_BEZIER_CURVE:
428  case VTK_BEZIER_TRIANGLE:
430  return 0;
431 
432  case VTK_TETRA:
433  case VTK_QUADRATIC_TETRA:
438  return 4;
439 
440  case VTK_PYRAMID:
444  case VTK_WEDGE:
445  case VTK_QUADRATIC_WEDGE:
449  case VTK_LAGRANGE_WEDGE:
450  case VTK_BEZIER_WEDGE:
451  return 5;
452 
453  case VTK_VOXEL:
454  case VTK_HEXAHEDRON:
462  return 6;
463 
465  return 7;
466 
467  case VTK_HEXAGONAL_PRISM:
468  return 8;
469 
470  case VTK_POLYHEDRON: // Need to look these up
471  if (!this->CheckCache(FacesFlag))
472  {
473  this->FetchFaces();
474  this->SetCache(FacesFlag);
475  }
476  return this->Faces->GetNumberOfIds() != 0 ? this->Faces->GetId(0) : 0;
477 
478  default:
479  vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
480  break;
481  }
482 
483  return 0;
484 }
485 
486 #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:140
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition: vtkIdList.h:166
vtkIdType GetId(const vtkIdType i)
Return the id at location i.
Definition: vtkIdList.h:171
a simple class to control print indentation
Definition: vtkIndent.h:113
abstract base class for most VTK objects
Definition: vtkObject.h:73
represent and manipulate 3D points
Definition: vtkPoints.h:143
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
@ VTK_VOXEL
Definition: vtkCellType.h:96
@ VTK_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:109
@ VTK_PARAMETRIC_SURFACE
Definition: vtkCellType.h:132
@ VTK_HIGHER_ORDER_TETRAHEDRON
Definition: vtkCellType.h:143
@ VTK_TRIANGLE_STRIP
Definition: vtkCellType.h:91
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:118
@ VTK_LAGRANGE_CURVE
Definition: vtkCellType.h:149
@ VTK_HIGHER_ORDER_QUAD
Definition: vtkCellType.h:141
@ VTK_PYRAMID
Definition: vtkCellType.h:99
@ VTK_PIXEL
Definition: vtkCellType.h:93
@ VTK_QUADRATIC_WEDGE
Definition: vtkCellType.h:110
@ VTK_BEZIER_WEDGE
Definition: vtkCellType.h:163
@ VTK_BIQUADRATIC_QUAD
Definition: vtkCellType.h:112
@ VTK_HIGHER_ORDER_WEDGE
Definition: vtkCellType.h:144
@ VTK_LAGRANGE_QUADRILATERAL
Definition: vtkCellType.h:151
@ VTK_POLY_LINE
Definition: vtkCellType.h:89
@ VTK_TRIQUADRATIC_PYRAMID
Definition: vtkCellType.h:114
@ VTK_TRIANGLE
Definition: vtkCellType.h:90
@ VTK_BEZIER_TRIANGLE
Definition: vtkCellType.h:159
@ VTK_POLYGON
Definition: vtkCellType.h:92
@ VTK_EMPTY_CELL
Definition: vtkCellType.h:85
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:111
@ VTK_POLYHEDRON
Definition: vtkCellType.h:128
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:113
@ VTK_TETRA
Definition: vtkCellType.h:95
@ VTK_LINE
Definition: vtkCellType.h:88
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:125
@ VTK_BEZIER_HEXAHEDRON
Definition: vtkCellType.h:162
@ VTK_PARAMETRIC_TRI_SURFACE
Definition: vtkCellType.h:133
@ VTK_LAGRANGE_WEDGE
Definition: vtkCellType.h:154
@ VTK_LAGRANGE_HEXAHEDRON
Definition: vtkCellType.h:153
@ VTK_PENTAGONAL_PRISM
Definition: vtkCellType.h:100
@ VTK_HIGHER_ORDER_TRIANGLE
Definition: vtkCellType.h:140
@ VTK_QUADRATIC_QUAD
Definition: vtkCellType.h:106
@ VTK_WEDGE
Definition: vtkCellType.h:98
@ VTK_PARAMETRIC_QUAD_SURFACE
Definition: vtkCellType.h:134
@ VTK_LAGRANGE_TETRAHEDRON
Definition: vtkCellType.h:152
@ VTK_PARAMETRIC_CURVE
Definition: vtkCellType.h:131
@ VTK_BEZIER_CURVE
Definition: vtkCellType.h:158
@ VTK_HIGHER_ORDER_PYRAMID
Definition: vtkCellType.h:145
@ VTK_HEXAGONAL_PRISM
Definition: vtkCellType.h:101
@ VTK_PARAMETRIC_HEX_REGION
Definition: vtkCellType.h:136
@ VTK_BEZIER_QUADRILATERAL
Definition: vtkCellType.h:160
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:116
@ VTK_HEXAHEDRON
Definition: vtkCellType.h:97
@ VTK_CUBIC_LINE
Definition: vtkCellType.h:122
@ VTK_LAGRANGE_TRIANGLE
Definition: vtkCellType.h:150
@ VTK_HIGHER_ORDER_HEXAHEDRON
Definition: vtkCellType.h:146
@ VTK_QUADRATIC_POLYGON
Definition: vtkCellType.h:107
@ VTK_QUAD
Definition: vtkCellType.h:94
@ VTK_QUADRATIC_TRIANGLE
Definition: vtkCellType.h:105
@ VTK_PARAMETRIC_TETRA_REGION
Definition: vtkCellType.h:135
@ VTK_QUADRATIC_EDGE
Definition: vtkCellType.h:104
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:108
@ VTK_HIGHER_ORDER_EDGE
Definition: vtkCellType.h:139
@ VTK_BEZIER_TETRAHEDRON
Definition: vtkCellType.h:161
@ VTK_VERTEX
Definition: vtkCellType.h:86
@ VTK_POLY_VERTEX
Definition: vtkCellType.h:87
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:115
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition: vtkCellType.h:117
@ VTK_HIGHER_ORDER_POLYGON
Definition: vtkCellType.h:142
@ VTK_BIQUADRATIC_TRIANGLE
Definition: vtkCellType.h:119
int vtkIdType
Definition: vtkType.h:332