VTK  9.2.6
vtkStaticPointLocator.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkStaticPointLocator.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
57#ifndef vtkStaticPointLocator_h
58#define vtkStaticPointLocator_h
59
61#include "vtkCommonDataModelModule.h" // For export macro
62
63class vtkIdList;
64struct vtkBucketList;
65class vtkDataArray;
66
67class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator : public vtkAbstractPointLocator
68{
69public:
75
77
81 void PrintSelf(ostream& os, vtkIndent indent) override;
83
85
90 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
91 vtkGetMacro(NumberOfPointsPerBucket, int);
93
95
101 vtkSetVector3Macro(Divisions, int);
102 vtkGetVectorMacro(Divisions, int, 3);
104
105 // Re-use any superclass signatures that we don't override.
110
118 vtkIdType FindClosestPoint(const double x[3]) override;
119
121
130 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
132 double radius, const double x[3], double inputDataLength, double& dist2);
134
143 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
144
151 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
152
162 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
163 double ptX[3], vtkIdType& ptId);
164
175 void MergePoints(double tol, vtkIdType* mergeMap);
176
188
190
194 void Initialize() override;
195 void FreeSearchStructure() override;
196 void BuildLocator() override;
197 void ForceBuildLocator() override;
198 void BuildLocator(const double* inBounds);
200
206 void GenerateRepresentation(int level, vtkPolyData* pd) override;
207
213
219 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
220
222
236 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
237 vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
239
247 bool GetLargeIds() { return this->LargeIds; }
248
250
254 virtual double* GetSpacing() { return this->H; }
255 virtual void GetSpacing(double spacing[3])
256 {
257 spacing[0] = this->H[0];
258 spacing[1] = this->H[1];
259 spacing[2] = this->H[2];
260 }
262
275 {
276 POINT_ORDER = 0,
277 BIN_ORDER = 1
278 };
279
281
286 vtkSetClampMacro(
288 vtkGetMacro(TraversalOrder, int);
290 {
291 this->SetTraversalOrder(vtkStaticPointLocator::POINT_ORDER);
292 }
295
296protected:
299
300 void BuildLocatorInternal() override;
301
302 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
303 int Divisions[3]; // Number of sub-divisions in x-y-z directions
304 double H[3]; // Width of each bucket in x-y-z directions
305 vtkBucketList* Buckets; // Lists of point ids in each bucket
306 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
307 bool LargeIds; // indicate whether integer ids are small or large
308 int TraversalOrder; // Control traversal order when threading
309
310private:
312 void operator=(const vtkStaticPointLocator&) = delete;
313};
314
315#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
list of point or cell ids
Definition vtkIdList.h:31
a simple class to control print indentation
Definition vtkIndent.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:85
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.
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.
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.
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 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.
int vtkIdType
Definition vtkType.h:332
#define VTK_ID_MAX
Definition vtkType.h:336
#define VTK_INT_MAX
Definition vtkType.h:155