VTK  9.6.1
vtkStaticPointLocator.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
58
59#ifndef vtkStaticPointLocator_h
60#define vtkStaticPointLocator_h
61
63#include "vtkCommonDataModelModule.h" // For export macro
64
65VTK_ABI_NAMESPACE_BEGIN
66class vtkDoubleArray;
67class vtkIdList;
68struct vtkBucketList;
69class vtkDataArray;
71
72class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator : public vtkAbstractPointLocator
73{
74public:
80
82
86 void PrintSelf(ostream& os, vtkIndent indent) override;
88
90
95 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
96 vtkGetMacro(NumberOfPointsPerBucket, int);
98
100
106 vtkSetVector3Macro(Divisions, int);
107 vtkGetVectorMacro(Divisions, int, 3);
109
110 // Reuse any superclass signatures that we don't override.
115
123 vtkIdType FindClosestPoint(const double x[3]) override;
124
126
135 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
137 double radius, const double x[3], double inputDataLength, double& dist2);
139
149 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
150
170 double FindNPointsInShell(int N, const double x[3], vtkDist2TupleArray& results,
171 double minDist2 = (-0.1), bool sort = true, vtkDoubleArray* petals = nullptr);
172
179 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
180
190 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
191 double ptX[3], vtkIdType& ptId);
192
203 void MergePoints(double tol, vtkIdType* mergeMap);
204
216
218
222 void Initialize() override;
223 void FreeSearchStructure() override;
224 void BuildLocator() override;
225 void ForceBuildLocator() override;
226 void BuildLocator(const double* inBounds);
228
234 void GenerateRepresentation(int level, vtkPolyData* pd) override;
235
241
247 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
248
253 void GetBucketCenter(int i, int j, int k, double center[3]);
254
256
270 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
273
281 bool GetLargeIds() { return this->LargeIds; }
282
284
288 virtual double* GetSpacing() { return this->H; }
289 virtual void GetSpacing(double spacing[3])
290 {
291 spacing[0] = this->H[0];
292 spacing[1] = this->H[1];
293 spacing[2] = this->H[2];
294 }
295
296
313
315
320 vtkSetClampMacro(
322 vtkGetMacro(TraversalOrder, int);
329
331
335 vtkSetMacro(Padding, double);
336 vtkGetMacro(Padding, double);
338
340
346 void StaticOn() { this->Static = true; }
347 void StaticOff() { this->Static = false; }
348 vtkGetMacro(Static, vtkTypeBool);
350
354 vtkBucketList* GetBuckets() { return this->Buckets; }
355
356protected:
359
360 void BuildLocatorInternal() override;
361
362 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
363 int Divisions[3]; // Number of sub-divisions in x-y-z directions
364 double H[3]; // Width of each bucket in x-y-z directions
365 vtkBucketList* Buckets; // Lists of point ids in each bucket
366 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
367 bool LargeIds; // Integer point ids are small (32-bit) or large (64-bit)
368 int TraversalOrder; // Control traversal order when threading
369 double Padding; // Pad out the bounding box of the locator
370 vtkTypeBool Static; // Control whether to repeatedly check modified time
371
372private:
374 void operator=(const vtkStaticPointLocator&) = delete;
375};
376
377VTK_ABI_NAMESPACE_END
378#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:24
a simple class to control print indentation
Definition vtkIndent.h:29
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:72
void StaticOff()
Turn on/off flag to control whether the locator checks modified time after it is built.
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.
virtual void SetTraversalOrder(int)
Specify the manner in which points are processed when a non-zero merge tolerance is specified.
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.
void GetBucketCenter(int i, int j, int k, double center[3])
Given a bucket/bin located at position (i,j,k), compute the center of the bucket.
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.
double FindNPointsInShell(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...
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 StaticOn()
Turn on/off flag to control whether the locator checks modified time after it is built.
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.
vtkBucketList * GetBuckets()
This method is useful for accessing the raw binned data.
Private declarations for 3D binned spatial locator.
Represent an array of vtkDist2Tuples.
int vtkTypeBool
Definition vtkABI.h:64
#define vtkDataArray
int vtkIdType
Definition vtkType.h:368
#define VTK_ID_MAX
Definition vtkType.h:372
#define VTK_INT_MAX
Definition vtkType.h:197