VTK  9.5.20251211
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
53
54#ifndef vtkStaticPointLocator2D_h
55#define vtkStaticPointLocator2D_h
56
58#include "vtkCommonDataModelModule.h" // For export macro
59
60VTK_ABI_NAMESPACE_BEGIN
61class vtkDoubleArray;
62class vtkIdList;
63struct vtkBucketList2D;
65
66class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator2D : public vtkAbstractPointLocator
67{
68public:
74
76
80 void PrintSelf(ostream& os, vtkIndent indent) override;
82
84
89 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
90 vtkGetMacro(NumberOfPointsPerBucket, int);
92
94
100 vtkSetVector2Macro(Divisions, int);
101 vtkGetVectorMacro(Divisions, int, 2);
103
104 // Reuse any superclass signatures that we don't override.
109
117 vtkIdType FindClosestPoint(const double x[3]) override;
118
120
128 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
130 double radius, const double x[3], double inputDataLength, double& dist2);
132
141 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
142
160 double FindNPointsInAnnulus(int N, const double x[3], vtkDist2TupleArray& results,
161 double minDist2 = (-0.1), bool sort = true, vtkDoubleArray* petals = nullptr);
162
169 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
170
180 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
181 double ptX[3], vtkIdType& ptId);
182
184
191 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
193
201 void MergePoints(double tol, vtkIdType* mergeMap);
202
204
208 void Initialize() override;
209 void FreeSearchStructure() override;
210 void BuildLocator() override;
211 void ForceBuildLocator() override;
213
221
227 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
228
230
244 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
247
255 bool GetLargeIds() { return this->LargeIds; }
256
258
261 void GetBounds(double* bounds) override
262 {
263 bounds[0] = this->Bounds[0];
264 bounds[1] = this->Bounds[1];
265 bounds[2] = this->Bounds[2];
266 bounds[3] = this->Bounds[3];
267 }
268
269
271
275 virtual double* GetSpacing() { return this->H; }
276 virtual void GetSpacing(double spacing[3])
277 {
278 spacing[0] = this->H[0];
279 spacing[1] = this->H[1];
280 spacing[2] = 0.0;
281 }
282
283
285
290 void GetBucketIndices(const double* x, int ij[2]) const;
291 vtkIdType GetBucketIndex(const double* x) const;
293
295
301 void StaticOn() { this->Static = true; }
302 void StaticOff() { this->Static = false; }
303 vtkGetMacro(Static, vtkTypeBool);
305
309 vtkBucketList2D* GetBuckets() { return this->Buckets; }
310
316 void GenerateRepresentation(int level, vtkPolyData* pd) override;
317
318protected:
321
322 void BuildLocatorInternal() override;
323
324 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
325 int Divisions[2]; // Number of sub-divisions in x-y directions
326 double H[2]; // Width of each bucket in x-y directions
327 vtkBucketList2D* Buckets; // Lists of point ids in each bucket
328 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
329 bool LargeIds; // indicate whether integer ids are small or large
330 vtkTypeBool Static; // Control whether to repeatedly check modified time
331
332private:
334 void operator=(const vtkStaticPointLocator2D&) = delete;
335};
336
337VTK_ABI_NAMESPACE_END
338#endif
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.
dynamic, self-adjusting array of double
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
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 StaticOn()
Turn on/off flag to control whether the locator checks modified time after it is built.
void GetBucketIndices(const double *x, int ij[2]) const
Given a point x[3], return the locator index (i,j) which contains the point.
double FindNPointsInAnnulus(int N, const double x[3], vtkDist2TupleArray &results, double minDist2=(-0.1), bool sort=true, vtkDoubleArray *petals=nullptr)
Find approximately N close points which are strictly greater than >minDist2 away from the query 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.
vtkBucketList2D * GetBuckets()
This method is useful for accessing the raw binned data.
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 StaticOff()
Turn on/off flag to control whether the locator checks modified time after it is built.
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...
Private declarations for 2D binned spatial locator.
Represent an array of vtkDist2Tuples.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:367
#define VTK_ID_MAX
Definition vtkType.h:371
#define VTK_INT_MAX
Definition vtkType.h:196