VTK
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkStaticPointLocator2D Class Reference

quickly locate points in 2-space More...

#include <vtkStaticPointLocator2D.h>

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

Public Member Functions

vtkIdType FindClosestPoint (const double x[3]) override
 Given a position x, return the id of the point closest to it. More...
 
void FindClosestNPoints (int N, const double x[3], vtkIdList *result) override
 Find the closest N points to a position. More...
 
void FindPointsWithinRadius (double R, const double x[3], vtkIdList *result) override
 Find all points within a specified radius R of position x. More...
 
int IntersectWithLine (double a0[3], double a1[3], double tol, double &t, double lineX[3], double ptX[3], vtkIdType &ptId)
 Intersect the points contained in the locator with the line defined by (a0,a1). More...
 
void MergePoints (double tol, vtkIdType *mergeMap)
 Merge points in the locator given a tolerance. More...
 
vtkIdType GetNumberOfPointsInBucket (vtkIdType bNum)
 Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return the number of points found in the bucket. More...
 
void GetBucketIds (vtkIdType bNum, vtkIdList *bList)
 Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return a list of point ids contained within the bucket. More...
 
bool GetLargeIds ()
 Inform the user as to whether large ids are being used. More...
 
void GenerateRepresentation (int level, vtkPolyData *pd) override
 Populate a polydata with the faces of the bins that potentially contain cells. More...
 
virtual void SetNumberOfPointsPerBucket (int)
 Specify the average number of points in each bucket. More...
 
virtual int GetNumberOfPointsPerBucket ()
 Specify the average number of points in each bucket. More...
 
virtual void SetDivisions (int, int)
 Set the number of divisions in x-y directions. More...
 
void SetDivisions (int[2])
 Set the number of divisions in x-y directions. More...
 
virtual intGetDivisions ()
 Set the number of divisions in x-y directions. More...
 
virtual void GetDivisions (int data[2])
 Set the number of divisions in x-y directions. More...
 
vtkIdType FindClosestPointWithinRadius (double radius, const double x[3], double &dist2) override
 Given a position x and a radius r, return the id of the point closest to the point within that radius. More...
 
virtual vtkIdType FindClosestPointWithinRadius (double radius, const double x[3], double inputDataLength, double &dist2)
 Given a position x and a radius r, return the id of the point closest to the point within that radius. More...
 
double FindCloseNBoundedPoints (int N, const double x[3], vtkIdList *result)
 Special method for 2D operations (e.g., vtkVoronoi2D). More...
 
void Initialize () override
 See vtkLocator and vtkAbstractPointLocator interface documentation. More...
 
void FreeSearchStructure () override
 See vtkLocator and vtkAbstractPointLocator interface documentation. More...
 
void BuildLocator () override
 See vtkLocator and vtkAbstractPointLocator interface documentation. More...
 
virtual void SetMaxNumberOfBuckets (vtkIdType)
 Set the maximum number of buckets in the locator. More...
 
virtual vtkIdType GetMaxNumberOfBuckets ()
 Set the maximum number of buckets in the locator. More...
 
void GetBounds (double *bounds) override
 Provide an accessor to the bounds. More...
 
virtual doubleGetSpacing ()
 Provide an accessor to the bucket spacing. More...
 
virtual void GetSpacing (double spacing[3])
 Provide an accessor to the bucket spacing. More...
 
void GetBucketIndices (const double *x, int ij[2]) const
 Given a point x[3], return the locator index (i,j) which contains the point. More...
 
vtkIdType GetBucketIndex (const double *x) const
 Given a point x[3], return the locator index (i,j) which contains the point. More...
 
- Public Member Functions inherited from vtkAbstractPointLocator
vtkIdType FindClosestPoint (double x, double y, double z)
 Given a position x, return the id of the point closest to it. More...
 
void FindClosestNPoints (int N, double x, double y, double z, vtkIdList *result)
 Find the closest N points to a position. More...
 
void FindPointsWithinRadius (double R, double x, double y, double z, vtkIdList *result)
 Find all points within a specified radius R of position x. More...
 
virtual doubleGetBounds ()
 Provide an accessor to the bounds. More...
 
virtual vtkIdType GetNumberOfBuckets ()
 Return the total number of buckets in the locator. More...
 
vtkAbstractPointLocatorNewInstance () const
 Standard type and print methods. More...
 
- Public Member Functions inherited from vtkLocator
virtual void Update ()
 Cause the locator to rebuild itself if it or its input dataset has changed. More...
 
virtual void SetDataSet (vtkDataSet *)
 Build the locator from the points/cells defining this dataset. More...
 
virtual vtkDataSetGetDataSet ()
 Build the locator from the points/cells defining this dataset. More...
 
virtual void SetMaxLevel (int)
 Set the maximum allowable level for the tree. More...
 
virtual int GetMaxLevel ()
 Set the maximum allowable level for the tree. More...
 
virtual int GetLevel ()
 Get the level of the locator (determined automatically if Automatic is true). More...
 
virtual void SetAutomatic (vtkTypeBool)
 Boolean controls whether locator depth/resolution of locator is computed automatically from average number of entities in bucket. More...
 
virtual vtkTypeBool GetAutomatic ()
 Boolean controls whether locator depth/resolution of locator is computed automatically from average number of entities in bucket. More...
 
virtual void AutomaticOn ()
 Boolean controls whether locator depth/resolution of locator is computed automatically from average number of entities in bucket. More...
 
virtual void AutomaticOff ()
 Boolean controls whether locator depth/resolution of locator is computed automatically from average number of entities in bucket. More...
 
virtual void SetTolerance (double)
 Specify absolute tolerance (in world coordinates) for performing geometric operations. More...
 
virtual double GetTolerance ()
 Specify absolute tolerance (in world coordinates) for performing geometric operations. More...
 
virtual vtkMTimeType GetBuildTime ()
 Return the time of the last data structure build. More...
 
void Register (vtkObjectBase *o) override
 Handle the PointSet <-> Locator loop. More...
 
void UnRegister (vtkObjectBase *o) override
 Handle the PointSet <-> Locator loop. More...
 
vtkLocatorNewInstance () const
 Standard type and print methods. More...
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on. More...
 
virtual void DebugOff ()
 Turn debugging output off. More...
 
bool GetDebug ()
 Get the value of the debug flag. More...
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag. More...
 
virtual void Modified ()
 Update the modification time for this object. More...
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time. More...
 
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)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
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. More...
 
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. More...
 
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. More...
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
virtual void Delete ()
 Delete a VTK object. More...
 
virtual void FastDelete ()
 Delete a reference to this object. More...
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream. More...
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
void PrintRevisions (ostream &)
 Legacy. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 

Static Public Member Functions

static vtkStaticPointLocator2DNew ()
 Construct with automatic computation of divisions, averaging 5 points per bucket. More...
 
- Static Public Member Functions inherited from vtkAbstractPointLocator
static vtkTypeBool IsTypeOf (const char *type)
 Standard type and print methods. More...
 
static vtkAbstractPointLocatorSafeDownCast (vtkObjectBase *o)
 Standard type and print methods. More...
 
- Static Public Member Functions inherited from vtkLocator
static vtkTypeBool IsTypeOf (const char *type)
 Standard type and print methods. More...
 
static vtkLocatorSafeDownCast (vtkObjectBase *o)
 Standard type and print methods. More...
 
- 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. More...
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes. More...
 
static void SetGlobalWarningDisplay (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static int GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
- 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. More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 

Protected Member Functions

 vtkStaticPointLocator2D ()
 
 ~vtkStaticPointLocator2D () override
 
- Protected Member Functions inherited from vtkAbstractPointLocator
 vtkAbstractPointLocator ()
 
 ~vtkAbstractPointLocator () override
 
- Protected Member Functions inherited from vtkLocator
 vtkLocator ()
 
 ~vtkLocator () override
 
void ReportReferences (vtkGarbageCollector *) override
 
- 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. More...
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events. More...
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

int NumberOfPointsPerBucket
 
int Divisions [2]
 
double H [2]
 
vtkBucketList2D * Buckets
 
vtkIdType MaxNumberOfBuckets
 
bool LargeIds
 
- Protected Attributes inherited from vtkAbstractPointLocator
double Bounds [6]
 
vtkIdType NumberOfBuckets
 
- Protected Attributes inherited from vtkLocator
vtkDataSetDataSet
 
vtkTypeBool Automatic
 
double Tolerance
 
int MaxLevel
 
int Level
 
vtkTimeStamp BuildTime
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 
typedef vtkAbstractPointLocator Superclass
 Standard type and print methods. More...
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard type and print methods. More...
 
static vtkStaticPointLocator2DSafeDownCast (vtkObjectBase *o)
 Standard type and print methods. More...
 
virtual vtkTypeBool IsA (const char *type)
 Standard type and print methods. More...
 
vtkStaticPointLocator2DNewInstance () const
 Standard type and print methods. More...
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard type and print methods. More...
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard type and print methods. More...
 

Additional Inherited Members

- Public Types inherited from vtkAbstractPointLocator
typedef vtkLocator Superclass
 Standard type and print methods. More...
 
- Public Types inherited from vtkLocator
typedef vtkObject Superclass
 Standard type and print methods. More...
 

Detailed Description

quickly locate points in 2-space

vtkStaticPointLocator2D is a spatial search object to quickly locate points in 2D. vtkStaticPointLocator2D works by dividing a specified region of space into a regular array of rectilinear buckets, and then keeping a list of points that lie in each bucket. Typical operation involves giving a position in 2D and finding the closest point; or finding the N closest points. (Note that the more general vtkStaticPointLocator is available for 3D operations.) Other specialized methods for 2D have also been provided.

vtkStaticPointLocator2D is an accelerated version of vtkPointLocator. It is threaded (via vtkSMPTools), and supports one-time static construction (i.e., incremental point insertion is not supported). If you need to incrementally insert points, use the vtkPointLocator or its kin to do so.

Note that to satisfy the superclass's API, methods often assume a 3D point is provided. However, only the x,y values are used for processing. The z-value is only used to define location of the 2D plane.

Warning
This class is templated. It may run slower than serial execution if the code is not optimized during compilation. Build in Release or ReleaseWithDebugInfo.
Make sure that you review the documentation for the superclasses vtkAbstactPointLocator and vtkLocator. In particular the Automatic data member can be used to automatically determine divisions based on the average number of points per bucket.
Other types of spatial locators have been developed such as octrees and kd-trees. These are often more efficient for the operations described here.
See also
vtkStaticPointLocator vtkPointLocator vtkCellLocator vtkLocator vtkAbstractPointLocator
Tests:
vtkStaticPointLocator2D (Tests)

Definition at line 69 of file vtkStaticPointLocator2D.h.

Member Typedef Documentation

Standard type and print methods.

Definition at line 82 of file vtkStaticPointLocator2D.h.

Constructor & Destructor Documentation

vtkStaticPointLocator2D::vtkStaticPointLocator2D ( )
protected
vtkStaticPointLocator2D::~vtkStaticPointLocator2D ( )
overrideprotected

Member Function Documentation

static vtkStaticPointLocator2D* vtkStaticPointLocator2D::New ( )
static

Construct with automatic computation of divisions, averaging 5 points per bucket.

static vtkTypeBool vtkStaticPointLocator2D::IsTypeOf ( const char *  type)
static

Standard type and print methods.

virtual vtkTypeBool vtkStaticPointLocator2D::IsA ( const char *  type)
virtual

Standard type and print methods.

Reimplemented from vtkAbstractPointLocator.

static vtkStaticPointLocator2D* vtkStaticPointLocator2D::SafeDownCast ( vtkObjectBase o)
static

Standard type and print methods.

virtual vtkObjectBase* vtkStaticPointLocator2D::NewInstanceInternal ( ) const
protectedvirtual

Standard type and print methods.

Reimplemented from vtkAbstractPointLocator.

vtkStaticPointLocator2D* vtkStaticPointLocator2D::NewInstance ( ) const

Standard type and print methods.

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

Standard type and print methods.

Reimplemented from vtkAbstractPointLocator.

virtual void vtkStaticPointLocator2D::SetNumberOfPointsPerBucket ( int  )
virtual

Specify the average number of points in each bucket.

This data member is used in conjunction with the Automatic data member (if enabled) to determine the number of locator x-y divisions.

virtual int vtkStaticPointLocator2D::GetNumberOfPointsPerBucket ( )
virtual

Specify the average number of points in each bucket.

This data member is used in conjunction with the Automatic data member (if enabled) to determine the number of locator x-y divisions.

virtual void vtkStaticPointLocator2D::SetDivisions ( int  ,
int   
)
virtual

Set the number of divisions in x-y directions.

If the Automatic data member is enabled, the Divisions are set according to the NumberOfPointsPerBucket and MaxNumberOfBuckets data members. The number of divisions must be >= 1 in each direction.

void vtkStaticPointLocator2D::SetDivisions ( int  [2])

Set the number of divisions in x-y directions.

If the Automatic data member is enabled, the Divisions are set according to the NumberOfPointsPerBucket and MaxNumberOfBuckets data members. The number of divisions must be >= 1 in each direction.

virtual int* vtkStaticPointLocator2D::GetDivisions ( )
virtual

Set the number of divisions in x-y directions.

If the Automatic data member is enabled, the Divisions are set according to the NumberOfPointsPerBucket and MaxNumberOfBuckets data members. The number of divisions must be >= 1 in each direction.

virtual void vtkStaticPointLocator2D::GetDivisions ( int  data[2])
virtual

Set the number of divisions in x-y directions.

If the Automatic data member is enabled, the Divisions are set according to the NumberOfPointsPerBucket and MaxNumberOfBuckets data members. The number of divisions must be >= 1 in each direction.

vtkIdType vtkStaticPointLocator2D::FindClosestPoint ( const double  x[3])
overridevirtual

Given a position x, return the id of the point closest to it.

An alternative method (defined in the superclass) requires separate x-y-z values. These methods are thread safe if BuildLocator() is directly or indirectly called from a single thread first. (Note in the 2D locator the z-value is ignored.)

Implements vtkAbstractPointLocator.

vtkIdType vtkStaticPointLocator2D::FindClosestPointWithinRadius ( double  radius,
const double  x[3],
double dist2 
)
overridevirtual

Given a position x and a radius r, return the id of the point closest to the point within that radius.

These methods are thread safe if BuildLocator() is directly or indirectly called from a single thread first. dist2 returns the squared distance to the point. Note that if multiple points are located the same distance away, the actual point returned is a function in which order the points are processed (i.e., indeterminate).

Implements vtkAbstractPointLocator.

virtual vtkIdType vtkStaticPointLocator2D::FindClosestPointWithinRadius ( double  radius,
const double  x[3],
double  inputDataLength,
double dist2 
)
virtual

Given a position x and a radius r, return the id of the point closest to the point within that radius.

These methods are thread safe if BuildLocator() is directly or indirectly called from a single thread first. dist2 returns the squared distance to the point. Note that if multiple points are located the same distance away, the actual point returned is a function in which order the points are processed (i.e., indeterminate).

void vtkStaticPointLocator2D::FindClosestNPoints ( int  N,
const double  x[3],
vtkIdList result 
)
overridevirtual

Find the closest N points to a position.

This returns the closest N points to a position. A faster method could be created that returned N close points to a position, but necessarily the exact N closest. The returned points are sorted from closest to farthest. These methods are thread safe if BuildLocator() is directly or indirectly called from a single thread first.

Implements vtkAbstractPointLocator.

void vtkStaticPointLocator2D::FindPointsWithinRadius ( double  R,
const double  x[3],
vtkIdList result 
)
overridevirtual

Find all points within a specified radius R of position x.

The result is not sorted in any specific manner. These methods are thread safe if BuildLocator() is directly or indirectly called from a single thread first.

Implements vtkAbstractPointLocator.

int vtkStaticPointLocator2D::IntersectWithLine ( double  a0[3],
double  a1[3],
double  tol,
double t,
double  lineX[3],
double  ptX[3],
vtkIdType ptId 
)

Intersect the points contained in the locator with the line defined by (a0,a1).

Return the point within the tolerance tol that is closest to a0 (tol measured in the world coordinate system). If an intersection occurs (i.e., the method returns nonzero), then the parametric location along the line t, the closest position along the line lineX, and the coordinates of the picked ptId is returned in ptX. (This method is thread safe after the locator is built.)

double vtkStaticPointLocator2D::FindCloseNBoundedPoints ( int  N,
const double  x[3],
vtkIdList result 
)

Special method for 2D operations (e.g., vtkVoronoi2D).

The method returns the approximate number of points requested, returning the radius R of the furthest point, with the guarantee that all points are included that are closer than <=R.

void vtkStaticPointLocator2D::MergePoints ( double  tol,
vtkIdType mergeMap 
)

Merge points in the locator given a tolerance.

Return a merge map which represents the mapping of "concident" point ids to a single point. Note the number of points in the merge map is the number of points the locator was built with. The user is expected to pass in an allocated mergeMap.

void vtkStaticPointLocator2D::Initialize ( )
overridevirtual

See vtkLocator and vtkAbstractPointLocator interface documentation.

These methods are not thread safe.

Reimplemented from vtkLocator.

void vtkStaticPointLocator2D::FreeSearchStructure ( )
overridevirtual

See vtkLocator and vtkAbstractPointLocator interface documentation.

These methods are not thread safe.

Implements vtkLocator.

void vtkStaticPointLocator2D::BuildLocator ( )
overridevirtual

See vtkLocator and vtkAbstractPointLocator interface documentation.

These methods are not thread safe.

Implements vtkLocator.

vtkIdType vtkStaticPointLocator2D::GetNumberOfPointsInBucket ( vtkIdType  bNum)

Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return the number of points found in the bucket.

Note that a bucket can also be specified with locator indices (i,j) which converts to a the bucket number bNum=i+this->Divisions[0]*j.

void vtkStaticPointLocator2D::GetBucketIds ( vtkIdType  bNum,
vtkIdList bList 
)

Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return a list of point ids contained within the bucket.

The user must provide an instance of vtkIdList to contain the result.

virtual void vtkStaticPointLocator2D::SetMaxNumberOfBuckets ( vtkIdType  )
virtual

Set the maximum number of buckets in the locator.

By default the value is set to VTK_INT_MAX. Note that there are significant performance implications at work here. If the number of buckets is set very large (meaning > VTK_INT_MAX) then internal sorting may be performed using 64-bit integers (which is much slower than using a 32-bit int). Of course, memory requirements may dramatically increase as well. It is recommended that the default value be used; but for extremely large data it may be desired to create a locator with an exceptionally large number of buckets. Note also that during initialization of the locator if the MaxNumberOfBuckets threshold is exceeded, the Divisions are scaled down in such a way as not to exceed the MaxNumberOfBuckets proportionally to the size of the bounding box in the x-y-z directions.

virtual vtkIdType vtkStaticPointLocator2D::GetMaxNumberOfBuckets ( )
virtual

Set the maximum number of buckets in the locator.

By default the value is set to VTK_INT_MAX. Note that there are significant performance implications at work here. If the number of buckets is set very large (meaning > VTK_INT_MAX) then internal sorting may be performed using 64-bit integers (which is much slower than using a 32-bit int). Of course, memory requirements may dramatically increase as well. It is recommended that the default value be used; but for extremely large data it may be desired to create a locator with an exceptionally large number of buckets. Note also that during initialization of the locator if the MaxNumberOfBuckets threshold is exceeded, the Divisions are scaled down in such a way as not to exceed the MaxNumberOfBuckets proportionally to the size of the bounding box in the x-y-z directions.

bool vtkStaticPointLocator2D::GetLargeIds ( )
inline

Inform the user as to whether large ids are being used.

This flag only has meaning after the locator has been built. Large ids are used when the number of binned points, or the number of bins, is >= the maximum number of buckets (specified by the user). Note that LargeIds are only available on 64-bit architectures.

Definition at line 240 of file vtkStaticPointLocator2D.h.

void vtkStaticPointLocator2D::GetBounds ( double bounds)
inlineoverridevirtual

Provide an accessor to the bounds.

Valid after the locator is built.

Reimplemented from vtkAbstractPointLocator.

Definition at line 246 of file vtkStaticPointLocator2D.h.

virtual double* vtkStaticPointLocator2D::GetSpacing ( )
inlinevirtual

Provide an accessor to the bucket spacing.

Valid after the locator is built.

Definition at line 260 of file vtkStaticPointLocator2D.h.

virtual void vtkStaticPointLocator2D::GetSpacing ( double  spacing[3])
inlinevirtual

Provide an accessor to the bucket spacing.

Valid after the locator is built.

Definition at line 261 of file vtkStaticPointLocator2D.h.

void vtkStaticPointLocator2D::GetBucketIndices ( const double x,
int  ij[2] 
) const

Given a point x[3], return the locator index (i,j) which contains the point.

This method is meant to be fast, so error checking is not performed. This method should only be called after the locator is built.

vtkIdType vtkStaticPointLocator2D::GetBucketIndex ( const double x) const

Given a point x[3], return the locator index (i,j) which contains the point.

This method is meant to be fast, so error checking is not performed. This method should only be called after the locator is built.

void vtkStaticPointLocator2D::GenerateRepresentation ( int  level,
vtkPolyData pd 
)
overridevirtual

Populate a polydata with the faces of the bins that potentially contain cells.

Note that the level parameter has no effect on this method as there is no hierarchy built (i.e., uniform binning). Typically this is used for debugging.

Implements vtkLocator.

Member Data Documentation

int vtkStaticPointLocator2D::NumberOfPointsPerBucket
protected

Definition at line 286 of file vtkStaticPointLocator2D.h.

int vtkStaticPointLocator2D::Divisions[2]
protected

Definition at line 287 of file vtkStaticPointLocator2D.h.

double vtkStaticPointLocator2D::H[2]
protected

Definition at line 288 of file vtkStaticPointLocator2D.h.

vtkBucketList2D* vtkStaticPointLocator2D::Buckets
protected

Definition at line 289 of file vtkStaticPointLocator2D.h.

vtkIdType vtkStaticPointLocator2D::MaxNumberOfBuckets
protected

Definition at line 290 of file vtkStaticPointLocator2D.h.

bool vtkStaticPointLocator2D::LargeIds
protected

Definition at line 291 of file vtkStaticPointLocator2D.h.


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