VTK  9.2.5
vtkOctreePointLocator.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkOctreePointLocator.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=========================================================================*/
15/*----------------------------------------------------------------------------
16 Copyright (c) Sandia Corporation
17 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18----------------------------------------------------------------------------*/
19
40#ifndef vtkOctreePointLocator_h
41#define vtkOctreePointLocator_h
42
44#include "vtkCommonDataModelModule.h" // For export macro
45
46class vtkCellArray;
47class vtkIdTypeArray;
49class vtkPoints;
50class vtkPolyData;
51
52class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
53{
54public:
56 void PrintSelf(ostream& os, vtkIndent indent) override;
57
59
61
64 vtkSetMacro(MaximumPointsPerRegion, int);
65 vtkGetMacro(MaximumPointsPerRegion, int);
67
69
72 vtkSetMacro(CreateCubicOctants, int);
73 vtkGetMacro(CreateCubicOctants, int);
75
77
83 vtkGetMacro(FudgeFactor, double);
84 vtkSetMacro(FudgeFactor, double);
86
88
92 double* GetBounds() override;
93 void GetBounds(double* bounds) override;
95
97
100 vtkGetMacro(NumberOfLeafNodes, int);
102
106 void GetRegionBounds(int regionID, double bounds[6]);
107
111 void GetRegionDataBounds(int leafNodeID, double bounds[6]);
112
116 int GetRegionContainingPoint(double x, double y, double z);
117
125 void BuildLocator() override;
126
130 void ForceBuildLocator() override;
131
133
137 vtkIdType FindClosestPoint(const double x[3]) override;
138 vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
140
146 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
147
149
154 vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
155 vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
157
162 void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
163
172 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
173
178
182 void FreeSearchStructure() override;
183
188 void GenerateRepresentation(int level, vtkPolyData* pd) override;
189
196 void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
197
198protected:
201
202 void BuildLocatorInternal() override;
203
205 vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
206
208
210
214 int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
215 int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
217
219
221
228 vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
229
230 // Recursive helper for public FindPointsWithinRadius
232
233 // Recursive helper for public FindPointsInArea
235
236 // Recursive helper for public FindPointsInArea
238
239 void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
240
241 int DivideTest(int size, int level);
242
244
249 int FindClosestPointInRegion_(int leafNodeId, double x, double y, double z, double& dist2);
250
259 double x, double y, double z, double radius, int skipRegion, double& dist2);
260
262
268
269 double FudgeFactor; // a very small distance, relative to the dataset's size
273
274 float MaxWidth;
275
284
286 void operator=(const vtkOctreePointLocator&) = delete;
287};
288#endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:187
list of point or cell ids
Definition: vtkIdList.h:34
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:40
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.
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).
static vtkOctreePointLocator * New()
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)
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.
vtkIdTypeArray * GetPointsInRegion(int leafNodeId)
Get a list of the original IDs of all points in a leaf node.
represent and manipulate 3D points
Definition: vtkPoints.h:40
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:91
int vtkIdType
Definition: vtkType.h:332
#define VTK_NEWINSTANCE