VTK  9.2.6
vtkStaticPointLocator.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkStaticPointLocator.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=========================================================================*/
60#ifndef vtkStaticPointLocator_h
61#define vtkStaticPointLocator_h
62
64#include "vtkCommonDataModelModule.h" // For export macro
65
66class vtkIdList;
67struct vtkBucketList;
68class vtkDataArray;
69
70class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator : public vtkAbstractPointLocator
71{
72public:
78
80
84 void PrintSelf(ostream& os, vtkIndent indent) override;
86
88
93 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
94 vtkGetMacro(NumberOfPointsPerBucket, int);
96
98
104 vtkSetVector3Macro(Divisions, int);
105 vtkGetVectorMacro(Divisions, int, 3);
107
108 // Re-use any superclass signatures that we don't override.
113
121 vtkIdType FindClosestPoint(const double x[3]) override;
122
124
133 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
135 double radius, const double x[3], double inputDataLength, double& dist2);
137
146 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
147
154 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
155
165 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
166 double ptX[3], vtkIdType& ptId);
167
178 void MergePoints(double tol, vtkIdType* mergeMap);
179
191
193
197 void Initialize() override;
198 void FreeSearchStructure() override;
199 void BuildLocator() override;
200 void ForceBuildLocator() override;
201 void BuildLocator(const double* inBounds);
203
209 void GenerateRepresentation(int level, vtkPolyData* pd) override;
210
216
222 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
223
225
239 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
240 vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
242
250 bool GetLargeIds() { return this->LargeIds; }
251
253
257 virtual double* GetSpacing() { return this->H; }
258 virtual void GetSpacing(double spacing[3])
259 {
260 spacing[0] = this->H[0];
261 spacing[1] = this->H[1];
262 spacing[2] = this->H[2];
263 }
265
278 {
279 POINT_ORDER = 0,
280 BIN_ORDER = 1
281 };
282
284
289 vtkSetClampMacro(
291 vtkGetMacro(TraversalOrder, int);
293 {
294 this->SetTraversalOrder(vtkStaticPointLocator::POINT_ORDER);
295 }
298
299protected:
302
303 void BuildLocatorInternal() override;
304
305 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
306 int Divisions[3]; // Number of sub-divisions in x-y-z directions
307 double H[3]; // Width of each bucket in x-y-z directions
308 vtkBucketList* Buckets; // Lists of point ids in each bucket
309 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
310 bool LargeIds; // indicate whether integer ids are small or large
311 int TraversalOrder; // Control traversal order when threading
312
313private:
315 void operator=(const vtkStaticPointLocator&) = delete;
316};
317
318#endif
abstract class to quickly locate points in 3-space
virtual void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)=0
Find all points within a specified radius R of position x.
virtual vtkIdType FindClosestPoint(const double x[3])=0
Given a position x, return the id of the point closest to it.
virtual double * GetBounds()
Provide an accessor to the bounds.
virtual void FindClosestNPoints(int N, const double x[3], vtkIdList *result)=0
Find the closest N points to a position.
abstract superclass for arrays of numeric data
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
quickly locate points in 3-space
void MergePoints(double tol, vtkIdType *mergeMap)
Merge points in the locator given a tolerance.
vtkIdType GetNumberOfPointsInBucket(vtkIdType bNum)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return the number of point...
void Initialize() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
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 in that radius,...
void SetTraversalOrderToBinOrder()
Specify the manner in which points are processed when a non-zero merge tolerance is specified.
void MergePointsWithData(vtkDataArray *data, vtkIdType *mergeMap)
Merge points and associated data in the locator.
void FindClosestNPoints(int N, const double x[3], vtkIdList *result) override
Find the closest N points to a position.
TraversalOrderType
Point merging is inherently an order-dependent process.
bool GetLargeIds()
Inform the user as to whether large ids are being used.
~vtkStaticPointLocator() override
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
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).
void FreeSearchStructure() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void SetTraversalOrderToPointOrder()
Specify the manner in which points are processed when a non-zero merge tolerance is specified.
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 in that radius,...
void ForceBuildLocator() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void GetBucketIds(vtkIdType bNum, vtkIdList *bList)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return a list of point ids...
void GenerateRepresentation(int level, vtkPolyData *pd) override
Populate a polydata with the faces of the bins that potentially contain cells.
void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result) override
Find all points within a specified radius R of position x.
static vtkStaticPointLocator * New()
Construct with automatic computation of divisions, averaging 5 points per bucket.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
vtkIdType FindClosestPoint(const double x[3]) override
Given a position x, return the id of the point closest to it, or (-1) if no point found.
void BuildLocator(const double *inBounds)
See vtkLocator and vtkAbstractPointLocator interface documentation.
void BuildLocator() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
virtual void GetSpacing(double spacing[3])
Provide an accessor to the bucket spacing.
int vtkIdType
Definition vtkType.h:332
#define VTK_ID_MAX
Definition vtkType.h:336
#define VTK_INT_MAX
Definition vtkType.h:155