VTK  9.3.0
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 
25 #ifndef vtkOctreePointLocator_h
26 #define vtkOctreePointLocator_h
27 
29 #include "vtkCommonDataModelModule.h" // For export macro
30 
31 VTK_ABI_NAMESPACE_BEGIN
32 class vtkCellArray;
33 class vtkIdTypeArray;
35 class vtkPoints;
36 class vtkPolyData;
37 
38 class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
39 {
40 public:
42  void PrintSelf(ostream& os, vtkIndent indent) override;
43 
45 
47 
50  vtkSetMacro(MaximumPointsPerRegion, int);
51  vtkGetMacro(MaximumPointsPerRegion, int);
53 
55 
58  vtkSetMacro(CreateCubicOctants, int);
59  vtkGetMacro(CreateCubicOctants, int);
61 
63 
69  vtkGetMacro(FudgeFactor, double);
70  vtkSetMacro(FudgeFactor, double);
72 
74 
78  double* GetBounds() override;
79  void GetBounds(double* bounds) override;
81 
83 
86  vtkGetMacro(NumberOfLeafNodes, int);
88 
92  void GetRegionBounds(int regionID, double bounds[6]);
93 
97  void GetRegionDataBounds(int leafNodeID, double bounds[6]);
98 
102  int GetRegionContainingPoint(double x, double y, double z);
103 
111  void BuildLocator() override;
112 
116  void ForceBuildLocator() override;
117 
119 
123  vtkIdType FindClosestPoint(const double x[3]) override;
124  vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
126 
132  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
133 
135 
140  vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
141  vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
143 
148  void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
149 
158  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
159 
164 
168  void FreeSearchStructure() override;
169 
174  void GenerateRepresentation(int level, vtkPolyData* pd) override;
175 
182  void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
183 
184 protected:
187 
188  void BuildLocatorInternal() override;
189 
191  vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
192 
194 
196 
200  int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
201  int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
203 
205 
207 
214  vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
215 
216  // Recursive helper for public FindPointsWithinRadius
218 
219  // Recursive helper for public FindPointsInArea
221 
222  // Recursive helper for public FindPointsInArea
224 
225  void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
226 
227  int DivideTest(int size, int level);
228 
230 
235  int FindClosestPointInRegion_(int leafNodeId, double x, double y, double z, double& dist2);
236 
245  double x, double y, double z, double radius, int skipRegion, double& dist2);
246 
248 
254 
255  double FudgeFactor; // a very small distance, relative to the dataset's size
259 
260  float MaxWidth;
261 
270 
272  void operator=(const vtkOctreePointLocator&) = delete;
273 };
274 VTK_ABI_NAMESPACE_END
275 #endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:176
list of point or cell ids
Definition: vtkIdList.h:23
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:29
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:29
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:80
@ 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