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
48 #ifndef vtkStaticPointLocator_h
49 #define vtkStaticPointLocator_h
50 
52 #include "vtkCommonDataModelModule.h" // For export macro
53 
54 VTK_ABI_NAMESPACE_BEGIN
55 class vtkIdList;
56 struct vtkBucketList;
57 class vtkDataArray;
58 
59 class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator : public vtkAbstractPointLocator
60 {
61 public:
67 
69 
73  void PrintSelf(ostream& os, vtkIndent indent) override;
75 
77 
82  vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
83  vtkGetMacro(NumberOfPointsPerBucket, int);
85 
87 
93  vtkSetVector3Macro(Divisions, int);
94  vtkGetVectorMacro(Divisions, int, 3);
96 
97  // Re-use any superclass signatures that we don't override.
102 
110  vtkIdType FindClosestPoint(const double x[3]) override;
111 
113 
122  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
124  double radius, const double x[3], double inputDataLength, double& dist2);
126 
135  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
136 
143  void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
144 
154  int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
155  double ptX[3], vtkIdType& ptId);
156 
167  void MergePoints(double tol, vtkIdType* mergeMap);
168 
180 
182 
186  void Initialize() override;
187  void FreeSearchStructure() override;
188  void BuildLocator() override;
189  void ForceBuildLocator() override;
190  void BuildLocator(const double* inBounds);
192 
198  void GenerateRepresentation(int level, vtkPolyData* pd) override;
199 
205 
211  void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
212 
214 
228  vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
229  vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
231 
239  bool GetLargeIds() { return this->LargeIds; }
240 
242 
246  virtual double* GetSpacing() { return this->H; }
247  virtual void GetSpacing(double spacing[3])
248  {
249  spacing[0] = this->H[0];
250  spacing[1] = this->H[1];
251  spacing[2] = this->H[2];
252  }
254 
267  {
268  POINT_ORDER = 0,
269  BIN_ORDER = 1
270  };
271 
273 
278  vtkSetClampMacro(
280  vtkGetMacro(TraversalOrder, int);
282  {
283  this->SetTraversalOrder(vtkStaticPointLocator::POINT_ORDER);
284  }
287 
288 protected:
291 
292  void BuildLocatorInternal() override;
293 
294  int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
295  int Divisions[3]; // Number of sub-divisions in x-y-z directions
296  double H[3]; // Width of each bucket in x-y-z directions
297  vtkBucketList* Buckets; // Lists of point ids in each bucket
298  vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
299  bool LargeIds; // indicate whether integer ids are small or large
300  int TraversalOrder; // Control traversal order when threading
301 
302 private:
304  void operator=(const vtkStaticPointLocator&) = delete;
305 };
306 
307 VTK_ABI_NAMESPACE_END
308 #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:45
list of point or cell ids
Definition: vtkIdList.h:23
a simple class to control print indentation
Definition: vtkIndent.h:29
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:80
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