VTK  9.4.20241014
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkSphericalPointIterator Class Reference

Traverse a collection of points in spherical ordering. More...

#include <vtkSphericalPointIterator.h>

Inheritance diagram for vtkSphericalPointIterator:
[legend]
Collaboration diagram for vtkSphericalPointIterator:
[legend]

Public Types

enum  AxesType {
  XY_CW_AXES = 0 , XY_CCW_AXES = 1 , XY_SQUARE_AXES = 2 , CUBE_AXES = 3 ,
  OCTAHEDRON_AXES = 4 , CUBE_OCTAHEDRON_AXES , DODECAHEDRON_AXES = 6 , ICOSAHEDRON_AXES = 7
}
 While the axes can be arbitrarily specified, it is possible to select axes from a menu of predefined axes sets. More...
 
enum  SortType { SORT_NONE = 0 , SORT_ASCENDING = 1 , SORT_DESCENDING = 2 }
 Points can be sorted along each axis. More...
 

Public Member Functions

void SetAxes (int axesType, int resolution=6)
 A convenience method to set the iterator axes from the predefined set enumerated above.
 
void GoToFirstPoint ()
 Begin iterating over the neighborhood of points.
 
bool IsDoneWithTraversal ()
 Return true if set traversal is completed.
 
void GoToNextPoint ()
 Go to the next point in the neighborhood.
 
void GetCurrentPoint (vtkIdType &ptId, double x[3])
 Get the current point (point id and coordinates) during forward iteration.
 
vtkIdType GetCurrentPoint ()
 Return the current point id during forward iteration.
 
vtkIdType GetPoint (int axis, int ptIdx)
 Provide random access to the jth point of the ith axis.
 
vtkIdType GetNumberOfAxes ()
 Return the number of axes defined.
 
void GetAxisPoints (int axis, vtkIdType &npts, const vtkIdType *&pts)
 Return the list of points along the specified ith axis.
 
void BuildRepresentation (vtkPolyData *pd)
 A convenience method that produces a geometric representation of the iterator (e.g., axes + center).
 
 vtkSetSmartPointerMacro (DataSet, vtkDataSet)
 Define the dataset and its associated points over which to iterate.
 
 vtkGetSmartPointerMacro (DataSet, vtkDataSet)
 Define the dataset and its associated points over which to iterate.
 
 vtkSetSmartPointerMacro (Axes, vtkDoubleArray)
 Define the axes for the point iterator.
 
 vtkGetSmartPointerMacro (Axes, vtkDoubleArray)
 Define the axes for the point iterator.
 
virtual void SetSorting (int)
 Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.
 
virtual int GetSorting ()
 Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.
 
void SetSortTypeToNone ()
 Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.
 
void SetSortTypeToAscending ()
 Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.
 
void SetSortTypeToDescending ()
 Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.
 
bool Initialize (double center[3], vtkIdList *neighborhood)
 Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.
 
bool Initialize (double center[3], vtkIdType numNei, vtkIdType *neighborhood)
 Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.
 
bool Initialize (double center[3])
 Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on.
 
virtual void DebugOff ()
 Turn debugging output off.
 
bool GetDebug ()
 Get the value of the debug flag.
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag.
 
virtual void Modified ()
 Update the modification time for this object.
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
vtkTypeBool HasObserver (unsigned long event)
 
vtkTypeBool HasObserver (const char *event)
 
vtkTypeBool InvokeEvent (unsigned long event)
 
vtkTypeBool InvokeEvent (const char *event)
 
std::string GetObjectDescription () const override
 The object description printed in messages and PrintSelf output.
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events.
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events.
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Allow user to set the AbortFlagOn() with the return value of the callback method.
 
vtkTypeBool InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
vtkTypeBool InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
virtual void SetObjectName (const std::string &objectName)
 Set/get the name of this object for reporting purposes.
 
virtual std::string GetObjectName () const
 Set/get the name of this object for reporting purposes.
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string.
 
virtual std::string GetObjectDescription () const
 The object description printed in messages and PrintSelf output.
 
virtual vtkTypeBool IsA (const char *name)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
virtual vtkIdType GetNumberOfGenerationsFromBase (const char *name)
 Given the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class).
 
virtual void Delete ()
 Delete a VTK object.
 
virtual void FastDelete ()
 Delete a reference to this object.
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream.
 
void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object).
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object).
 
int GetReferenceCount ()
 Return the current reference count of this object.
 
void SetReferenceCount (int)
 Sets the reference count.
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space.
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 
virtual bool UsesGarbageCollector () const
 Indicate whether the class uses vtkGarbageCollector or not.
 

Protected Member Functions

 vtkSphericalPointIterator ()
 
 ~vtkSphericalPointIterator () override=default
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
 ~vtkObject () override
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr)
 These methods allow a command to exclusively grab all events.
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events.
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void RegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
virtual void ObjectFinalize ()
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkSmartPointer< vtkDataSetDataSet
 
vtkSmartPointer< vtkDoubleArrayAxes
 
int Sorting
 
std::unique_ptr< SphericalPointIterator > Iterator
 
vtkTimeStamp BuildTime
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
std::string ObjectName
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 
static vtkSphericalPointIteratorNew ()
 Standard methods to instantiate, obtain type information, and print information about an instance of the class.
 
 vtkAbstractTypeMacro (vtkSphericalPointIterator, vtkObject)
 Standard methods to instantiate, obtain type information, and print information about an instance of the class.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods to instantiate, obtain type information, and print information about an instance of the class.
 

Additional Inherited Members

- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes.
 
static void SetGlobalWarningDisplay (vtkTypeBool val)
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static vtkTypeBool GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 Return 1 if this class type is the same type of (or a subclass of) the named class.
 
static vtkIdType GetNumberOfGenerationsFromBaseType (const char *name)
 Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class).
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within.
 
static bool GetUsingMemkind ()
 A global state flag that controls whether vtkObjects are constructed in the usual way (the default) or within the extended memory space.
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 

Detailed Description

Traverse a collection of points in spherical ordering.

vtkSphericalPointIterator is a state-based iterator for traversing a set of points (i.e., a neighborhood of points) in a dataset, providing a point traversal order across user-defined "axes" which span a 2D or 3D space (typically a circle or sphere). The points along each axes may be sorted in increasing radial order. To define the points, specify a dataset (i.e., its associated points, whether the points are represented implicitly or explicitly) and an associated neighborhood over which to iterate. Methods for iterating over the points are provided.

For example, consider the axes of iteration to be the four rays emanating from the center of a square and passing through the center of each of the four edges of the square. Points to be iterated over are associated (using a dot product) with each of the four axes, and then can be sorted along each axis. Then the order of iteration is then: (axis0,pt0), (axis1,pt0), (axis2,pt0), (axis3,pt0), (axis0,pt1), (axis1,pt1), (axis2,pt1), (axis3,pt1), (axis0,pt2), (axis1,pt2), (axis2,pt2), (axis3,pt2), and so on in a "spiraling" fashion until all points are visited. Thus the order of visitation is: iteration i visits all N axes in order, returning the jth point sorted along each of the N axes (i.e., i increases the fastest). Alternatively, methods exist to randomly access points, or points associated with an axes, so that custom iteration methods can be defined.

The iterator can be defined with any number of axes (defined by 3D vectors). The axes must not be coincident, and typically are equally spaced from one another. The order which the axes are defined determines the order in which the axes (and hence the points) are traversed. So for example, in a 2D sphere, four axes in the (-x,+x,-y,+y) directions would provide a "ping pong" iteration, while four axes ordered in the (+x,+y,-x,-y) directions would provide a counterclockwise rotation iteration.

The iterator provides thread-safe iteration of dataset points. It supports both random and forward iteration.

Warning
The behavior of the iterator depends on the ordering of the iteration axes. It is possible to obtain a wide variety of iteration patterns depending on these axes. For example, if only one axis is defined, then a "linear" pattern is possible (i.e., visiting points in the half space defined by the vector); if two axes, then a "diagonal" iteration pattern; and so on. Note that points are sorted along the iteration axes depending on the their projection onto them (e.g., using the dot product). Because only points with positive projection are associated with an axis, it is possible that some points in the neighborhood will not be processed (i.e., if a point in the neighborhood does not positively project onto any of the axes, then it will not be iterated over). Thus if all points are to be iterated over, then the axes must form a basis which covers all points using positive projections.
See also
vtkVoronoi2D vtkVoronoi3D vtkStaticPointLocator vtkPointLocator
Tests:
vtkSphericalPointIterator (Tests)

Definition at line 81 of file vtkSphericalPointIterator.h.

Member Enumeration Documentation

◆ AxesType

While the axes can be arbitrarily specified, it is possible to select axes from a menu of predefined axes sets.

For example, XY_CW_AXES refer to axes that rotate in clockwise direction starting with the first axis parallel to the x-axis; with the total number of axes given by the resolution. Some 3D regular polyhedral solids are also referred to: the axes pass through the center of the faces of the solid. So DODECAHEDRON axes refer to the 12 axes that pass through the center of the 12 faces of the dodecahedron. In some cases the resolution parameter need not be specified.

Enumerator
XY_CW_AXES 
XY_CCW_AXES 
XY_SQUARE_AXES 
CUBE_AXES 
OCTAHEDRON_AXES 
CUBE_OCTAHEDRON_AXES 
DODECAHEDRON_AXES 
ICOSAHEDRON_AXES 

Definition at line 131 of file vtkSphericalPointIterator.h.

◆ SortType

Points can be sorted along each axis.

By default, no sorting is performed. Other options are ascending and descending radial order. Ascended sorting results in point traversal starting near the center of the iterator, and proceeding radially outward. Descended sorting results in point traversal starting away from the center of the iterator, and proceeding radially inward.

Enumerator
SORT_NONE 
SORT_ASCENDING 
SORT_DESCENDING 

Definition at line 160 of file vtkSphericalPointIterator.h.

Constructor & Destructor Documentation

◆ vtkSphericalPointIterator()

vtkSphericalPointIterator::vtkSphericalPointIterator ( )
protected

◆ ~vtkSphericalPointIterator()

vtkSphericalPointIterator::~vtkSphericalPointIterator ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static vtkSphericalPointIterator * vtkSphericalPointIterator::New ( )
static

Standard methods to instantiate, obtain type information, and print information about an instance of the class.

◆ vtkAbstractTypeMacro()

vtkSphericalPointIterator::vtkAbstractTypeMacro ( vtkSphericalPointIterator  ,
vtkObject   
)

Standard methods to instantiate, obtain type information, and print information about an instance of the class.

◆ PrintSelf()

void vtkSphericalPointIterator::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

Standard methods to instantiate, obtain type information, and print information about an instance of the class.

Reimplemented from vtkObject.

◆ vtkSetSmartPointerMacro() [1/2]

vtkSphericalPointIterator::vtkSetSmartPointerMacro ( DataSet  ,
vtkDataSet   
)

Define the dataset and its associated points over which to iterate.

◆ vtkGetSmartPointerMacro() [1/2]

vtkSphericalPointIterator::vtkGetSmartPointerMacro ( DataSet  ,
vtkDataSet   
)

Define the dataset and its associated points over which to iterate.

◆ vtkSetSmartPointerMacro() [2/2]

vtkSphericalPointIterator::vtkSetSmartPointerMacro ( Axes  ,
vtkDoubleArray   
)

Define the axes for the point iterator.

This only needs to be defined once (typically immediately after instantiation). The axes data array must be a 3-component array, where each 3-tuple defines a vector defining a axis. The number of axes is limited to 100,000 or less (typically the numbers of axes are <=20). The order in which the axes are defined determines the order in which the axes are traversed. Depending on the order, it's possible to create a variety of traversal patterns, ranging from clockwise/counterclockwise to spiral to ping pong (e.g., -x,+x, -y,+y, -z,+z). Note: the defining axes need not be normalized, they are normalized and copied into internal iterator storage in the Initialize() method.

◆ vtkGetSmartPointerMacro() [2/2]

vtkSphericalPointIterator::vtkGetSmartPointerMacro ( Axes  ,
vtkDoubleArray   
)

Define the axes for the point iterator.

This only needs to be defined once (typically immediately after instantiation). The axes data array must be a 3-component array, where each 3-tuple defines a vector defining a axis. The number of axes is limited to 100,000 or less (typically the numbers of axes are <=20). The order in which the axes are defined determines the order in which the axes are traversed. Depending on the order, it's possible to create a variety of traversal patterns, ranging from clockwise/counterclockwise to spiral to ping pong (e.g., -x,+x, -y,+y, -z,+z). Note: the defining axes need not be normalized, they are normalized and copied into internal iterator storage in the Initialize() method.

◆ SetAxes()

void vtkSphericalPointIterator::SetAxes ( int  axesType,
int  resolution = 6 
)

A convenience method to set the iterator axes from the predefined set enumerated above.

The resolution parameter is optional in some cases - it is used by axes types that are non-fixed such as rotation of a vector around a center point in the plane (e.g., x-y plane).

◆ SetSorting()

virtual void vtkSphericalPointIterator::SetSorting ( int  )
virtual

Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.

(Note that some operators such as the locator query FindClosestNPoints() return radially sorted neighborhoods in ascending direction and often do not need sorting - this can save significant time.)

◆ GetSorting()

virtual int vtkSphericalPointIterator::GetSorting ( )
virtual

Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.

(Note that some operators such as the locator query FindClosestNPoints() return radially sorted neighborhoods in ascending direction and often do not need sorting - this can save significant time.)

◆ SetSortTypeToNone()

void vtkSphericalPointIterator::SetSortTypeToNone ( )
inline

Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.

(Note that some operators such as the locator query FindClosestNPoints() return radially sorted neighborhoods in ascending direction and often do not need sorting - this can save significant time.)

Definition at line 177 of file vtkSphericalPointIterator.h.

◆ SetSortTypeToAscending()

void vtkSphericalPointIterator::SetSortTypeToAscending ( )
inline

Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.

(Note that some operators such as the locator query FindClosestNPoints() return radially sorted neighborhoods in ascending direction and often do not need sorting - this can save significant time.)

Definition at line 178 of file vtkSphericalPointIterator.h.

◆ SetSortTypeToDescending()

void vtkSphericalPointIterator::SetSortTypeToDescending ( )
inline

Specify whether points along each axis are radially sorted, and if so, whether in an ascending or descending direction.

(Note that some operators such as the locator query FindClosestNPoints() return radially sorted neighborhoods in ascending direction and often do not need sorting - this can save significant time.)

Definition at line 179 of file vtkSphericalPointIterator.h.

◆ Initialize() [1/3]

bool vtkSphericalPointIterator::Initialize ( double  center[3],
vtkIdList neighborhood 
)

Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.

(The point ids refer to the points contained in the dataset.) If initialization fails (because the Axes or the DataSet have not been defined) then false is returned; true otherwise. One of the Initialize() variants enables iteration over all points in the dataset.

◆ Initialize() [2/3]

bool vtkSphericalPointIterator::Initialize ( double  center[3],
vtkIdType  numNei,
vtkIdType neighborhood 
)

Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.

(The point ids refer to the points contained in the dataset.) If initialization fails (because the Axes or the DataSet have not been defined) then false is returned; true otherwise. One of the Initialize() variants enables iteration over all points in the dataset.

◆ Initialize() [3/3]

bool vtkSphericalPointIterator::Initialize ( double  center[3])

Initialize the iteration process around a position [x], over a set of points (the neighborhood) defined by a list of numNei point ids.

(The point ids refer to the points contained in the dataset.) If initialization fails (because the Axes or the DataSet have not been defined) then false is returned; true otherwise. One of the Initialize() variants enables iteration over all points in the dataset.

◆ GoToFirstPoint()

void vtkSphericalPointIterator::GoToFirstPoint ( )

Begin iterating over the neighborhood of points.

It is possible that not all points are iterated over - those points not projecting onto any axis with a positive dot product are not visited.

◆ IsDoneWithTraversal()

bool vtkSphericalPointIterator::IsDoneWithTraversal ( )

Return true if set traversal is completed.

Otherwise false.

◆ GoToNextPoint()

void vtkSphericalPointIterator::GoToNextPoint ( )

Go to the next point in the neighborhood.

This is only valid when IsDoneWithTraversal() returns false;

◆ GetCurrentPoint() [1/2]

void vtkSphericalPointIterator::GetCurrentPoint ( vtkIdType ptId,
double  x[3] 
)

Get the current point (point id and coordinates) during forward iteration.

◆ GetCurrentPoint() [2/2]

vtkIdType vtkSphericalPointIterator::GetCurrentPoint ( )

Return the current point id during forward iteration.

◆ GetPoint()

vtkIdType vtkSphericalPointIterator::GetPoint ( int  axis,
int  ptIdx 
)

Provide random access to the jth point of the ith axis.

Returns the point id located at (axis,ptIdx); or a value <0 if the requested point does not exist.

◆ GetNumberOfAxes()

vtkIdType vtkSphericalPointIterator::GetNumberOfAxes ( )

Return the number of axes defined.

The value returned is valid only after Initialize() is invoked.

◆ GetAxisPoints()

void vtkSphericalPointIterator::GetAxisPoints ( int  axis,
vtkIdType npts,
const vtkIdType *&  pts 
)

Return the list of points along the specified ith axis.

◆ BuildRepresentation()

void vtkSphericalPointIterator::BuildRepresentation ( vtkPolyData pd)

A convenience method that produces a geometric representation of the iterator (e.g., axes + center).

The representation simply draws lines for each of the axes emanating from the center point. Each line (or line cell) is assigned cell data which is the axis number. This is typically used for debugging or educational purposes. Note that the method is valid only after Initialize() has been invoked.

Member Data Documentation

◆ DataSet

vtkSmartPointer<vtkDataSet> vtkSphericalPointIterator::DataSet
protected

Definition at line 261 of file vtkSphericalPointIterator.h.

◆ Axes

vtkSmartPointer<vtkDoubleArray> vtkSphericalPointIterator::Axes
protected

Definition at line 262 of file vtkSphericalPointIterator.h.

◆ Sorting

int vtkSphericalPointIterator::Sorting
protected

Definition at line 263 of file vtkSphericalPointIterator.h.

◆ Iterator

std::unique_ptr<SphericalPointIterator> vtkSphericalPointIterator::Iterator
protected

Definition at line 267 of file vtkSphericalPointIterator.h.

◆ BuildTime

vtkTimeStamp vtkSphericalPointIterator::BuildTime
protected

Definition at line 270 of file vtkSphericalPointIterator.h.


The documentation for this class was generated from the following file: