VTK  9.3.0
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
57 #ifndef vtkStaticPointLocator_h
58 #define vtkStaticPointLocator_h
59 
61 #include "vtkCommonDataModelModule.h" // For export macro
62 
63 VTK_ABI_NAMESPACE_BEGIN
64 class vtkIdList;
65 struct vtkBucketList;
66 class vtkDataArray;
67 
68 class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator : public vtkAbstractPointLocator
69 {
70 public:
76 
78 
82  void PrintSelf(ostream& os, vtkIndent indent) override;
84 
86 
91  vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
92  vtkGetMacro(NumberOfPointsPerBucket, int);
94 
96 
102  vtkSetVector3Macro(Divisions, int);
103  vtkGetVectorMacro(Divisions, int, 3);
105 
106  // Re-use any superclass signatures that we don't override.
111 
119  vtkIdType FindClosestPoint(const double x[3]) override;
120 
122 
131  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
133  double radius, const double x[3], double inputDataLength, double& dist2);
135 
144  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
145 
152  void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
153 
163  int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
164  double ptX[3], vtkIdType& ptId);
165 
176  void MergePoints(double tol, vtkIdType* mergeMap);
177 
189 
191 
195  void Initialize() override;
196  void FreeSearchStructure() override;
197  void BuildLocator() override;
198  void ForceBuildLocator() override;
199  void BuildLocator(const double* inBounds);
201 
207  void GenerateRepresentation(int level, vtkPolyData* pd) override;
208 
214 
220  void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
221 
223 
237  vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
238  vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
240 
248  bool GetLargeIds() { return this->LargeIds; }
249 
251 
255  virtual double* GetSpacing() { return this->H; }
256  virtual void GetSpacing(double spacing[3])
257  {
258  spacing[0] = this->H[0];
259  spacing[1] = this->H[1];
260  spacing[2] = this->H[2];
261  }
263 
276  {
277  POINT_ORDER = 0,
278  BIN_ORDER = 1
279  };
280 
282 
287  vtkSetClampMacro(
289  vtkGetMacro(TraversalOrder, int);
291  {
292  this->SetTraversalOrder(vtkStaticPointLocator::POINT_ORDER);
293  }
296 
297 protected:
300 
301  void BuildLocatorInternal() override;
302 
303  int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
304  int Divisions[3]; // Number of sub-divisions in x-y-z directions
305  double H[3]; // Width of each bucket in x-y-z directions
306  vtkBucketList* Buckets; // Lists of point ids in each bucket
307  vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
308  bool LargeIds; // indicate whether integer ids are small or large
309  int TraversalOrder; // Control traversal order when threading
310 
311 private:
313  void operator=(const vtkStaticPointLocator&) = delete;
314 };
315 
316 VTK_ABI_NAMESPACE_END
317 #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
Definition: vtkDataArray.h:54
list of point or cell ids
Definition: vtkIdList.h:32
a simple class to control print indentation
Definition: vtkIndent.h:38
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:89
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.
static vtkStaticPointLocator * New()
Construct with automatic computation of divisions, averaging 5 points per bucket.
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.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
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.
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
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.
@ level
Definition: vtkX3D.h:395
@ spacing
Definition: vtkX3D.h:481
@ radius
Definition: vtkX3D.h:252
@ data
Definition: vtkX3D.h:315
int vtkIdType
Definition: vtkType.h:315
#define VTK_ID_MAX
Definition: vtkType.h:319
#define VTK_INT_MAX
Definition: vtkType.h:144