VTK  9.2.20230320
vtkSphericalPointIterator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkSphericalPointIterator.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 
77 #ifndef vtkSphericalPointIterator_h
78 #define vtkSphericalPointIterator_h
79 
80 #include "vtkCommonDataModelModule.h" // For export macro
81 #include "vtkDataSet.h" // the dataset and its points to iterate over
82 #include "vtkDoubleArray.h" // For axes
83 #include "vtkObject.h"
84 #include "vtkSmartPointer.h" // auto destruct
85 
86 #include <memory> // for std::unique_ptr
87 
88 VTK_ABI_NAMESPACE_BEGIN
89 class vtkDoubleArray;
90 class vtkPolyData;
91 struct SpiralPointIterator;
92 
93 class VTKCOMMONDATAMODEL_EXPORT vtkSphericalPointIterator : public vtkObject
94 {
95 public:
97 
103  void PrintSelf(ostream& os, vtkIndent indent) override;
105 
107 
113 
115 
131 
143  enum AxesType
144  {
145  XY_CW_AXES = 0, // axes clockwise around center in x-y plane (resolution required)
146  XY_CCW_AXES = 1, // axes counterclockwise around center (resolution required)
147  XY_SQUARE_AXES = 2, // axes +x,-x, +y,-y: axes through the four faces of a square
148  CUBE_AXES = 3, // axes +x,-x, +y,-y, +z,-z: axes through the six faces of a cube
149  OCTAHEDRON_AXES = 4, // axes through the eight faces of a regular octahedron
150  CUBE_OCTAHEDRON_AXES =
151  5, // axes through the eight faces of a regular octahedron and six faces of a cube
152  DODECAHEDRON_AXES = 6, // axes through the twelve faces of a dedecahdron
153  ICOSAHEDRON_AXES = 7, // axes through the twenty faces of a icosahedron
154  };
155 
162  void SetAxes(int axesType, int resolution = 6);
163 
172  enum SortType
173  {
174  SORT_NONE = 0,
175  SORT_ASCENDING = 1,
176  SORT_DESCENDING = 2
177  };
178 
180 
187  vtkSetClampMacro(Sorting, int, SORT_NONE, SORT_DESCENDING);
188  vtkGetMacro(Sorting, int);
189  void SetSortTypeToNone() { this->SetSorting(SORT_NONE); }
190  void SetSortTypeToAscending() { this->SetSorting(SORT_ASCENDING); }
191  void SetSortTypeToDescending() { this->SetSorting(SORT_DESCENDING); }
193 
194  // The following methods support point iteration. The data members referred
195  // to previously must be defined before these iteration methods can be
196  // successfully invoked.
197 
199 
207  bool Initialize(double center[3], vtkIdList* neighborhood);
208  bool Initialize(double center[3], vtkIdType numNei, vtkIdType* neighborhood);
209  bool Initialize(double center[3]); // all points of the specified dataset
211 
218 
223 
229 
234  void GetCurrentPoint(vtkIdType& ptId, double x[3]);
235 
240 
245  vtkIdType GetPoint(int axis, int ptIdx);
246 
252 
256  void GetAxisPoints(int axis, vtkIdType& npts, const vtkIdType*& pts) VTK_SIZEHINT(pts, npts);
257 
267 
268 protected:
270  ~vtkSphericalPointIterator() override = default;
271 
272  // Information needed to define the spherical iterator.
273  vtkSmartPointer<vtkDataSet> DataSet; // The points to iterate over
274  vtkSmartPointer<vtkDoubleArray> Axes; // The axes defining the iteration pattern
275  int Sorting; // The direction of sorting, if sorting required
276 
277  // Iterator internals are represented using a PIMPL idiom
278  struct SphericalPointIterator;
279  std::unique_ptr<SphericalPointIterator> Iterator;
280 
281  // Changes to the VTK class must be propagated to the internal iterator
283 
284 private:
286  void operator=(const vtkSphericalPointIterator&) = delete;
287 };
288 
289 VTK_ABI_NAMESPACE_END
290 #endif // vtkSphericalPointIterator_h
abstract class to specify dataset behavior
Definition: vtkDataSet.h:174
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:144
a simple class to control print indentation
Definition: vtkIndent.h:120
abstract base class for most VTK objects
Definition: vtkObject.h:83
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:201
Traverse a collection of points in spherical ordering.
bool IsDoneWithTraversal()
Return true if set traversal is completed.
void SetSortTypeToNone()
Specify whether points along each axis are radially sorted, and if so, whether in an ascending or des...
SortType
Points can be sorted along each axis.
void SetAxes(int axesType, int resolution=6)
A convenience method to set the iterator axes from the predefined set enumerated above.
vtkSmartPointer< vtkDoubleArray > Axes
void SetSortTypeToDescending()
Specify whether points along each axis are radially sorted, and if so, whether in an ascending or des...
void GetCurrentPoint(vtkIdType &ptId, double x[3])
Get the current point (point id and coordinates) during forward iteration.
static vtkSphericalPointIterator * New()
Standard methods to instantiate, obtain type information, and print information about an instance of ...
vtkSetSmartPointerMacro(Axes, vtkDoubleArray)
Define the axes for the point iterator.
void GoToFirstPoint()
Begin iterating over the neighborhood of points.
void GetAxisPoints(int axis, vtkIdType &npts, const vtkIdType *&pts)
Return the list of points along the specified ith axis.
vtkGetSmartPointerMacro(DataSet, vtkDataSet)
Define the dataset and its associated points over which to iterate.
vtkSmartPointer< vtkDataSet > DataSet
std::unique_ptr< SphericalPointIterator > Iterator
bool Initialize(double center[3])
Initialize the iteration process around a position [x], over a set of points (the neighborhood) defin...
AxesType
While the axes can be arbitrarily specified, it is possible to select axes from a menu of predefined ...
void SetSortTypeToAscending()
Specify whether points along each axis are radially sorted, and if so, whether in an ascending or des...
void BuildRepresentation(vtkPolyData *pd)
A convenience method that produces a geometric representation of the iterator (e.g....
void GoToNextPoint()
Go to the next point in the neighborhood.
vtkSetSmartPointerMacro(DataSet, vtkDataSet)
Define the dataset and its associated points over which to iterate.
~vtkSphericalPointIterator() override=default
vtkIdType GetPoint(int axis, int ptIdx)
Provide random access to the jth point of the ith axis.
vtkGetSmartPointerMacro(Axes, vtkDoubleArray)
Define the axes for the point iterator.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to instantiate, obtain type information, and print information about an instance of ...
vtkAbstractTypeMacro(vtkSphericalPointIterator, vtkObject)
Standard methods to instantiate, obtain type information, and print information about an instance of ...
bool Initialize(double center[3], vtkIdType numNei, vtkIdType *neighborhood)
Initialize the iteration process around a position [x], over a set of points (the neighborhood) defin...
bool Initialize(double center[3], vtkIdList *neighborhood)
Initialize the iteration process around a position [x], over a set of points (the neighborhood) defin...
vtkIdType GetNumberOfAxes()
Return the number of axes defined.
vtkIdType GetCurrentPoint()
Return the current point id during forward iteration.
record modification and/or execution time
Definition: vtkTimeStamp.h:56
@ resolution
Definition: vtkX3D.h:478
@ center
Definition: vtkX3D.h:242
std::map< std::string, DataArray > DataSet
key: variable name, value: DataArray
Definition: VTXTypes.h:40
int vtkIdType
Definition: vtkType.h:327
#define VTK_SIZEHINT(...)