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
73
74#ifndef vtkSphericalPointIterator_h
75#define vtkSphericalPointIterator_h
76
77#include "vtkCommonDataModelModule.h" // For export macro
78#include "vtkDataSet.h" // the dataset and its points to iterate over
79#include "vtkDoubleArray.h" // For axes
80#include "vtkObject.h"
81#include "vtkSmartPointer.h" // auto destruct
82
83#include <memory> // for std::unique_ptr
84
85class vtkDoubleArray;
86class vtkPolyData;
87struct SpiralPointIterator;
88
89class VTKCOMMONDATAMODEL_EXPORT vtkSphericalPointIterator : public vtkObject
90{
91public:
93
99 void PrintSelf(ostream& os, vtkIndent indent) override;
101
103
109
111
127
140 {
141 XY_CW_AXES = 0, // axes clockwise around center in x-y plane (resolution required)
142 XY_CCW_AXES = 1, // axes counterclockwise around center (resolution required)
143 XY_SQUARE_AXES = 2, // axes +x,-x, +y,-y: axes through the four faces of a square
144 CUBE_AXES = 3, // axes +x,-x, +y,-y, +z,-z: axes through the six faces of a cube
145 OCTAHEDRON_AXES = 4, // axes through the eight faces of a regular octahedron
147 5, // axes through the eight faces of a regular octahedron and six faces of a cube
148 DODECAHEDRON_AXES = 6, // axes through the twelve faces of a dedecahdron
149 ICOSAHEDRON_AXES = 7, // axes through the twenty faces of a icosahedron
150 };
151
158 void SetAxes(int axesType, int resolution = 6);
159
174
176
183 vtkSetClampMacro(Sorting, int, SORT_NONE, SORT_DESCENDING);
184 vtkGetMacro(Sorting, int);
189
190 // The following methods support point iteration. The data members referred
191 // to previously must be defined before these iteration methods can be
192 // successfully invoked.
193
195
203 bool Initialize(double center[3], vtkIdList* neighborhood);
204 bool Initialize(double center[3], vtkIdType numNei, vtkIdType* neighborhood);
205 bool Initialize(double center[3]); // all points of the specified dataset
207
214
219
225
230 void GetCurrentPoint(vtkIdType& ptId, double x[3]);
231
236
241 vtkIdType GetPoint(int axis, int ptIdx);
242
248
252 void GetAxisPoints(int axis, vtkIdType& npts, const vtkIdType*& pts) VTK_SIZEHINT(pts, npts);
253
263
264protected:
266 ~vtkSphericalPointIterator() override = default;
267
268 // Information needed to define the spherical iterator.
269 vtkSmartPointer<vtkDataSet> DataSet; // The points to iterate over
270 vtkSmartPointer<vtkDoubleArray> Axes; // The axes defining the iteration pattern
271 int Sorting; // The direction of sorting, if sorting required
272
273 // Iterator internals are represented using a PIMPL idiom
274 struct SphericalPointIterator;
275 std::unique_ptr<SphericalPointIterator> Iterator;
276
277 // Changes to the VTK class must be propagated to the internal iterator
279
280private:
282 void operator=(const vtkSphericalPointIterator&) = delete;
283};
284
285#endif // vtkSphericalPointIterator_h
abstract class to specify dataset behavior
Definition vtkDataSet.h:57
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:31
a simple class to control print indentation
Definition vtkIndent.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:85
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(...)