VTK  9.2.5
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 
167 #ifndef vtkCellIterator_h
168 #define vtkCellIterator_h
169 
170 #include "vtkCellType.h" // For VTK_EMPTY_CELL
171 #include "vtkCommonDataModelModule.h" // For export macro
172 #include "vtkIdList.h" // For inline methods
173 #include "vtkNew.h" // For vtkNew
174 #include "vtkObject.h"
175 
176 class vtkGenericCell;
177 class vtkPoints;
178 
179 class VTKCOMMONDATAMODEL_EXPORT vtkCellIterator : public vtkObject
180 {
181 public:
182  void PrintSelf(ostream& os, vtkIndent indent) override;
184 
188  void InitTraversal();
189 
193  void GoToNextCell();
194 
198  virtual bool IsDoneWithTraversal() = 0;
199 
204  int GetCellType();
205 
211 
215  virtual vtkIdType GetCellId() = 0;
216 
221  vtkIdList* GetPointIds();
222 
228  vtkPoints* GetPoints();
229 
234  vtkIdList* GetFaces();
235 
241  void GetCell(vtkGenericCell* cell);
242 
247  vtkIdType GetNumberOfPoints();
248 
253  vtkIdType GetNumberOfFaces();
254 
255 protected:
257  ~vtkCellIterator() override;
258 
262  virtual void ResetToFirstCell() = 0;
263 
267  virtual void IncrementToNextCell() = 0;
268 
272  virtual void FetchCellType() = 0;
273 
277  virtual void FetchPointIds() = 0;
278 
282  virtual void FetchPoints() = 0;
283 
290  virtual void FetchFaces() {}
291 
292  int CellType;
296 
297 private:
298  vtkCellIterator(const vtkCellIterator&) = delete;
299  void operator=(const vtkCellIterator&) = delete;
300 
301  enum
302  {
303  UninitializedFlag = 0x0,
304  CellTypeFlag = 0x1,
305  PointIdsFlag = 0x2,
306  PointsFlag = 0x4,
307  FacesFlag = 0x8
308  };
309 
310  void ResetCache()
311  {
312  this->CacheFlags = UninitializedFlag;
313  this->CellType = VTK_EMPTY_CELL;
314  }
315 
316  void SetCache(unsigned char flags) { this->CacheFlags |= flags; }
317 
318  bool CheckCache(unsigned char flags) { return (this->CacheFlags & flags) == flags; }
319 
320  vtkNew<vtkPoints> PointsContainer;
321  vtkNew<vtkIdList> PointIdsContainer;
322  vtkNew<vtkIdList> FacesContainer;
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 //------------------------------------------------------------------------------
386 {
387  if (!this->CheckCache(PointIdsFlag))
388  {
389  this->FetchPointIds();
390  this->SetCache(PointIdsFlag);
391  }
392  return this->PointIds->GetNumberOfIds();
393 }
394 
395 //------------------------------------------------------------------------------
397 {
398  switch (this->GetCellType())
399  {
400  case VTK_EMPTY_CELL:
401  case VTK_VERTEX:
402  case VTK_POLY_VERTEX:
403  case VTK_LINE:
404  case VTK_POLY_LINE:
405  case VTK_TRIANGLE:
406  case VTK_TRIANGLE_STRIP:
407  case VTK_POLYGON:
408  case VTK_PIXEL:
409  case VTK_QUAD:
410  case VTK_QUADRATIC_EDGE:
412  case VTK_QUADRATIC_QUAD:
417  case VTK_CUBIC_LINE:
427  case VTK_LAGRANGE_CURVE:
430  case VTK_BEZIER_CURVE:
431  case VTK_BEZIER_TRIANGLE:
433  return 0;
434 
435  case VTK_TETRA:
436  case VTK_QUADRATIC_TETRA:
441  return 4;
442 
443  case VTK_PYRAMID:
447  case VTK_WEDGE:
448  case VTK_QUADRATIC_WEDGE:
452  case VTK_LAGRANGE_WEDGE:
453  case VTK_BEZIER_WEDGE:
454  return 5;
455 
456  case VTK_VOXEL:
457  case VTK_HEXAHEDRON:
465  return 6;
466 
468  return 7;
469 
470  case VTK_HEXAGONAL_PRISM:
471  return 8;
472 
473  case VTK_POLYHEDRON: // Need to look these up
474  if (!this->CheckCache(FacesFlag))
475  {
476  this->FetchFaces();
477  this->SetCache(FacesFlag);
478  }
479  return this->Faces->GetNumberOfIds() != 0 ? this->Faces->GetId(0) : 0;
480 
481  default:
482  vtkGenericWarningMacro("Unknown cell type: " << this->CellType);
483  break;
484  }
485 
486  return 0;
487 }
488 
489 #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:143
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition: vtkIdList.h:169
vtkIdType GetId(const vtkIdType i)
Return the id at location i.
Definition: vtkIdList.h:174
a simple class to control print indentation
Definition: vtkIndent.h:119
abstract base class for most VTK objects
Definition: vtkObject.h:82
represent and manipulate 3D points
Definition: vtkPoints.h:149
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