VTK  9.2.6
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
76
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
88class vtkDoubleArray;
89class vtkPolyData;
90struct SpiralPointIterator;
91
92class VTKCOMMONDATAMODEL_EXPORT vtkSphericalPointIterator : public vtkObject
93{
94public:
96
102 void PrintSelf(ostream& os, vtkIndent indent) override;
104
106
112
114
130
143 {
144 XY_CW_AXES = 0, // axes clockwise around center in x-y plane (resolution required)
145 XY_CCW_AXES = 1, // axes counterclockwise around center (resolution required)
146 XY_SQUARE_AXES = 2, // axes +x,-x, +y,-y: axes through the four faces of a square
147 CUBE_AXES = 3, // axes +x,-x, +y,-y, +z,-z: axes through the six faces of a cube
148 OCTAHEDRON_AXES = 4, // axes through the eight faces of a regular octahedron
150 5, // axes through the eight faces of a regular octahedron and six faces of a cube
151 DODECAHEDRON_AXES = 6, // axes through the twelve faces of a dedecahdron
152 ICOSAHEDRON_AXES = 7, // axes through the twenty faces of a icosahedron
153 };
154
161 void SetAxes(int axesType, int resolution = 6);
162
177
179
186 vtkSetClampMacro(Sorting, int, SORT_NONE, SORT_DESCENDING);
187 vtkGetMacro(Sorting, int);
192
193 // The following methods support point iteration. The data members referred
194 // to previously must be defined before these iteration methods can be
195 // successfully invoked.
196
198
206 bool Initialize(double center[3], vtkIdList* neighborhood);
207 bool Initialize(double center[3], vtkIdType numNei, vtkIdType* neighborhood);
208 bool Initialize(double center[3]); // all points of the specified dataset
210
217
222
228
233 void GetCurrentPoint(vtkIdType& ptId, double x[3]);
234
239
244 vtkIdType GetPoint(int axis, int ptIdx);
245
251
255 void GetAxisPoints(int axis, vtkIdType& npts, const vtkIdType*& pts) VTK_SIZEHINT(pts, npts);
256
266
267protected:
269 ~vtkSphericalPointIterator() override = default;
270
271 // Information needed to define the spherical iterator.
272 vtkSmartPointer<vtkDataSet> DataSet; // The points to iterate over
273 vtkSmartPointer<vtkDoubleArray> Axes; // The axes defining the iteration pattern
274 int Sorting; // The direction of sorting, if sorting required
275
276 // Iterator internals are represented using a PIMPL idiom
277 struct SphericalPointIterator;
278 std::unique_ptr<SphericalPointIterator> Iterator;
279
280 // Changes to the VTK class must be propagated to the internal iterator
282
283private:
285 void operator=(const vtkSphericalPointIterator&) = delete;
286};
287
288#endif // vtkSphericalPointIterator_h
abstract class to specify dataset behavior
Definition vtkDataSet.h:63
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:34
a simple class to control print indentation
Definition vtkIndent.h:40
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:91
Hold a reference to a vtkObjectBase instance.
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.
virtual void SetSorting(int)
Specify whether points along each axis are radially sorted, and if so, whether in an ascending or des...
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...
static vtkSphericalPointIterator * New()
Standard methods to instantiate, obtain type information, and print information about an instance of ...
vtkIdType GetNumberOfAxes()
Return the number of axes defined.
vtkIdType GetCurrentPoint()
Return the current point id during forward iteration.
record modification and/or execution time
int vtkIdType
Definition vtkType.h:332
#define VTK_SIZEHINT(...)