VTK  9.3.20240506
vtkStaticPointLocator2D.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
47#ifndef vtkStaticPointLocator2D_h
48#define vtkStaticPointLocator2D_h
49
51#include "vtkCommonDataModelModule.h" // For export macro
52
53VTK_ABI_NAMESPACE_BEGIN
54class vtkIdList;
55struct vtkBucketList2D;
56
57class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator2D : public vtkAbstractPointLocator
58{
59public:
65
67
71 void PrintSelf(ostream& os, vtkIndent indent) override;
73
75
80 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
81 vtkGetMacro(NumberOfPointsPerBucket, int);
83
85
91 vtkSetVector2Macro(Divisions, int);
92 vtkGetVectorMacro(Divisions, int, 2);
94
95 // Re-use any superclass signatures that we don't override.
100
108 vtkIdType FindClosestPoint(const double x[3]) override;
109
111
119 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
121 double radius, const double x[3], double inputDataLength, double& dist2);
123
132 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
133
140 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
141
151 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
152 double ptX[3], vtkIdType& ptId);
153
155
161 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
163
171 void MergePoints(double tol, vtkIdType* mergeMap);
172
174
178 void Initialize() override;
179 void FreeSearchStructure() override;
180 void BuildLocator() override;
181 void ForceBuildLocator() override;
183
191
197 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
198
200
214 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
215 vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
217
225 bool GetLargeIds() { return this->LargeIds; }
226
228
231 void GetBounds(double* bounds) override
232 {
233 bounds[0] = this->Bounds[0];
234 bounds[1] = this->Bounds[1];
235 bounds[2] = this->Bounds[2];
236 bounds[3] = this->Bounds[3];
237 }
239
241
245 virtual double* GetSpacing() { return this->H; }
246 virtual void GetSpacing(double spacing[3])
247 {
248 spacing[0] = this->H[0];
249 spacing[1] = this->H[1];
250 spacing[2] = 0.0;
251 }
253
255
260 void GetBucketIndices(const double* x, int ij[2]) const;
261 vtkIdType GetBucketIndex(const double* x) const;
263
269 void GenerateRepresentation(int level, vtkPolyData* pd) override;
270
271protected:
274
275 void BuildLocatorInternal() override;
276
277 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
278 int Divisions[2]; // Number of sub-divisions in x-y directions
279 double H[2]; // Width of each bucket in x-y directions
280 vtkBucketList2D* Buckets; // Lists of point ids in each bucket
281 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
282 bool LargeIds; // indicate whether integer ids are small or large
283
284private:
286 void operator=(const vtkStaticPointLocator2D&) = delete;
287};
288
289VTK_ABI_NAMESPACE_END
290#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.
list of point or cell ids
Definition vtkIdList.h:133
a simple class to control print indentation
Definition vtkIndent.h:108
concrete dataset represents vertices, lines, polygons, and triangle strips
quickly locate points in 2-space
vtkIdType FindClosestPoint(const double x[3]) override
Given a position x, return the id of the point closest to it.
void BuildLocator() 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 within that radius...
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
void FindClosestNPoints(int N, const double x[3], vtkIdList *result) override
Find the closest N points to a position.
double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList *result)
Special method for 2D operations (e.g., vtkVoronoi2D).
~vtkStaticPointLocator2D() override
bool GetLargeIds()
Inform the user as to whether large ids are being used.
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
void GetBucketIndices(const double *x, int ij[2]) const
Given a point x[3], return the locator index (i,j) which contains the point.
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...
vtkIdType GetBucketIndex(const double *x) const
Given a point x[3], return the locator index (i,j) which contains the point.
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 Initialize() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result) override
Find all points within a specified radius R of position x.
void GetBounds(double *bounds) override
Provide an accessor to the bounds.
void GenerateRepresentation(int level, vtkPolyData *pd) override
Populate a polydata with the faces of the bins that potentially contain cells.
void GetBucketIds(vtkIdType bNum, vtkIdList *bList)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return a list of point ids...
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
static vtkStaticPointLocator2D * New()
Construct with automatic computation of divisions, averaging 5 points per bucket.
void ForceBuildLocator() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void FreeSearchStructure() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void MergePoints(double tol, vtkIdType *mergeMap)
Merge points in the locator given a tolerance.
virtual void GetSpacing(double spacing[3])
Provide an accessor to the bucket spacing.
vtkIdType GetNumberOfPointsInBucket(vtkIdType bNum)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return the number of point...
int vtkIdType
Definition vtkType.h:315
#define VTK_ID_MAX
Definition vtkType.h:319
#define VTK_INT_MAX
Definition vtkType.h:144