VTK  9.3.20240327
vtkOctreePointLocator.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3 // SPDX-License-Identifier: BSD-3-Clause
4 
114 #ifndef vtkOctreePointLocator_h
115 #define vtkOctreePointLocator_h
116 
117 #include "vtkAbstractPointLocator.h"
118 #include "vtkCommonDataModelModule.h" // For export macro
119 
120 VTK_ABI_NAMESPACE_BEGIN
121 class vtkCellArray;
122 class vtkIdTypeArray;
124 class vtkPoints;
125 class vtkPolyData;
126 
127 class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
128 {
129 public:
131  void PrintSelf(ostream& os, vtkIndent indent) override;
132 
134 
136 
139  vtkSetMacro(MaximumPointsPerRegion, int);
140  vtkGetMacro(MaximumPointsPerRegion, int);
142 
144 
147  vtkSetMacro(CreateCubicOctants, int);
148  vtkGetMacro(CreateCubicOctants, int);
150 
152 
158  vtkGetMacro(FudgeFactor, double);
159  vtkSetMacro(FudgeFactor, double);
161 
163 
167  double* GetBounds() override;
168  void GetBounds(double* bounds) override;
170 
172 
175  vtkGetMacro(NumberOfLeafNodes, int);
177 
181  void GetRegionBounds(int regionID, double bounds[6]);
182 
186  void GetRegionDataBounds(int leafNodeID, double bounds[6]);
187 
191  int GetRegionContainingPoint(double x, double y, double z);
192 
200  void BuildLocator() override;
201 
205  void ForceBuildLocator() override;
206 
208 
212  vtkIdType FindClosestPoint(const double x[3]) override;
213  vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
215 
221  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
222 
224 
229  vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
230  vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
232 
237  void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
238 
247  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
248 
253 
257  void FreeSearchStructure() override;
258 
263  void GenerateRepresentation(int level, vtkPolyData* pd) override;
264 
271  void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
272 
273 protected:
276 
277  void BuildLocatorInternal() override;
278 
280  vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
281 
283 
285 
289  int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
290  int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
292 
294 
296 
303  vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
304 
305  // Recursive helper for public FindPointsWithinRadius
307 
308  // Recursive helper for public FindPointsInArea
310 
311  // Recursive helper for public FindPointsInArea
313 
314  void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
315 
316  int DivideTest(int size, int level);
317 
319 
324  int FindClosestPointInRegion_(int leafNodeId, double x, double y, double z, double& dist2);
325 
334  double x, double y, double z, double radius, int skipRegion, double& dist2);
335 
337 
343 
344  double FudgeFactor; // a very small distance, relative to the dataset's size
348 
349  float MaxWidth;
350 
359 
361  void operator=(const vtkOctreePointLocator&) = delete;
362 };
363 VTK_ABI_NAMESPACE_END
364 #endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:285
list of point or cell ids
Definition: vtkIdList.h:132
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:108
Octree node that has 8 children each of equal size.
an octree spatial decomposition of a set of points
void AddPolys(vtkOctreePointLocatorNode *node, vtkPoints *pts, vtkCellArray *polys)
void FindPointsWithinRadius(double radius, const double x[3], vtkIdList *result) override
Find all points within a specified radius of position x.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
int FindClosestPointInSphere(double x, double y, double z, double radius, int skipRegion, double &dist2)
Given a location and a radiues, find the closest point within this radius.
vtkIdTypeArray * GetPointsInRegion(int leafNodeId)
Get a list of the original IDs of all points in a leaf node.
vtkIdType FindClosestPointInRegion(int regionId, double *x, double &dist2)
Find the Id of the point in the given leaf region which is closest to the given point.
void ForceBuildLocator() override
Build the locator from the input dataset (even if UseExistingSearchStructure is on).
void FindClosestNPoints(int N, const double x[3], vtkIdList *result) override
Find the closest N points to a position.
void FindPointsInArea(vtkOctreePointLocatorNode *node, double *area, vtkIdTypeArray *ids)
void AddAllPointsInRegion(vtkOctreePointLocatorNode *node, vtkIdTypeArray *ids)
static vtkOctreePointLocator * New()
vtkOctreePointLocator(const vtkOctreePointLocator &)=delete
static void DeleteAllDescendants(vtkOctreePointLocatorNode *octant)
int DivideTest(int size, int level)
int NumberOfLeafNodes
The maximum number of points in a region/octant before it is subdivided.
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 FindPointsInArea(double *area, vtkIdTypeArray *ids, bool clearArray=true)
Fill ids with points found in area.
double * GetBounds() override
Get the spatial bounds of the entire octree space.
void BuildLeafNodeList(vtkOctreePointLocatorNode *node, int &index)
void AddAllPointsInRegion(vtkOctreePointLocatorNode *node, vtkIdList *ids)
vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double &dist2)
Find the Id of the point in the given leaf region which is closest to the given point.
int FindRegion(vtkOctreePointLocatorNode *node, double x, double y, double z)
Given a point and a node return the leaf node id that contains the point.
void GetRegionDataBounds(int leafNodeID, double bounds[6])
Get the bounds of the data within the leaf node.
int FindRegion(vtkOctreePointLocatorNode *node, float x, float y, float z)
Given a point and a node return the leaf node id that contains the point.
static void SetDataBoundsToSpatialBounds(vtkOctreePointLocatorNode *node)
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
vtkIdType FindClosestPoint(double x, double y, double z, double &dist2)
Return the Id of the point that is closest to the given point.
int FindClosestPointInRegion_(int leafNodeId, double x, double y, double z, double &dist2)
Given a leaf node id and point, return the local id and the squared distance between the closest poin...
vtkOctreePointLocatorNode * Top
~vtkOctreePointLocator() override
vtkIdType FindClosestPoint(const double x[3]) override
Return the Id of the point that is closest to the given point.
void GenerateRepresentation(int level, vtkPolyData *pd) override
Create a polydata representation of the boundaries of the octree regions.
int GetRegionContainingPoint(double x, double y, double z)
Get the id of the leaf region containing the specified location.
int MaximumPointsPerRegion
The maximum number of points in a region/octant before it is subdivided.
void DivideRegion(vtkOctreePointLocatorNode *node, int *ordering, int level)
int CreateCubicOctants
If CreateCubicOctants is non-zero, the bounding box of the points will be expanded such that all octa...
vtkOctreePointLocatorNode ** LeafNodeList
void GetRegionBounds(int regionID, double bounds[6])
Get the spatial bounds of octree region.
void BuildLocator() override
Create the octree decomposition of the cells of the data set or data sets.
void operator=(const vtkOctreePointLocator &)=delete
void GetBounds(double *bounds) override
Get the spatial bounds of the entire octree space.
void FreeSearchStructure() override
Delete the octree data structure.
void FindPointsWithinRadius(vtkOctreePointLocatorNode *node, double radiusSquared, const double x[3], vtkIdList *ids)
Recursive helper for public FindPointsWithinRadius.
represent and manipulate 3D points
Definition: vtkPoints.h:138
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:180
@ level
Definition: vtkX3D.h:395
@ radius
Definition: vtkX3D.h:252
@ size
Definition: vtkX3D.h:253
@ index
Definition: vtkX3D.h:246
int vtkIdType
Definition: vtkType.h:315
#define VTK_NEWINSTANCE