VTK  9.1.0
vtkStaticPointLocator2D.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkStaticPointLocator2D.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=========================================================================*/
56#ifndef vtkStaticPointLocator2D_h
57#define vtkStaticPointLocator2D_h
58
60#include "vtkCommonDataModelModule.h" // For export macro
61
62class vtkIdList;
63struct vtkBucketList2D;
64
65class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator2D : public vtkAbstractPointLocator
66{
67public:
73
75
79 void PrintSelf(ostream& os, vtkIndent indent) override;
81
83
88 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
89 vtkGetMacro(NumberOfPointsPerBucket, int);
91
93
99 vtkSetVector2Macro(Divisions, int);
100 vtkGetVectorMacro(Divisions, int, 2);
102
103 // Re-use any superclass signatures that we don't override.
108
116 vtkIdType FindClosestPoint(const double x[3]) override;
117
119
127 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
129 double radius, const double x[3], double inputDataLength, double& dist2);
131
140 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
141
148 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
149
159 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
160 double ptX[3], vtkIdType& ptId);
161
163
169 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
171
179 void MergePoints(double tol, vtkIdType* mergeMap);
180
182
186 void Initialize() override;
187 void FreeSearchStructure() override;
188 void BuildLocator() override;
190
198
204 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
205
207
221 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
222 vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
224
232 bool GetLargeIds() { return this->LargeIds; }
233
235
238 void GetBounds(double* bounds) override
239 {
240 bounds[0] = this->Bounds[0];
241 bounds[1] = this->Bounds[1];
242 bounds[2] = this->Bounds[2];
243 bounds[3] = this->Bounds[3];
244 }
246
248
252 virtual double* GetSpacing() { return this->H; }
253 virtual void GetSpacing(double spacing[3])
254 {
255 spacing[0] = this->H[0];
256 spacing[1] = this->H[1];
257 spacing[2] = 0.0;
258 }
260
262
267 void GetBucketIndices(const double* x, int ij[2]) const;
268 vtkIdType GetBucketIndex(const double* x) const;
270
276 void GenerateRepresentation(int level, vtkPolyData* pd) override;
277
278protected:
281
282 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
283 int Divisions[2]; // Number of sub-divisions in x-y directions
284 double H[2]; // Width of each bucket in x-y directions
285 vtkBucketList2D* Buckets; // Lists of point ids in each bucket
286 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
287 bool LargeIds; // indicate whether integer ids are small or large
288
289private:
291 void operator=(const vtkStaticPointLocator2D&) = delete;
292};
293
294#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.
list of point or cell ids
Definition: vtkIdList.h:140
a simple class to control print indentation
Definition: vtkIndent.h:113
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
quickly locate points in 2-space
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 GetBucketIndices(const double *x, int ij[2]) const
Given a point x[3], return the locator index (i,j) which contains the 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.
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 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...
@ level
Definition: vtkX3D.h:401
@ spacing
Definition: vtkX3D.h:487
@ radius
Definition: vtkX3D.h:258
int vtkIdType
Definition: vtkType.h:332
#define VTK_ID_MAX
Definition: vtkType.h:336
#define VTK_INT_MAX
Definition: vtkType.h:155