VTK  9.1.0
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
46#ifndef vtkOctreePointLocator_h
47#define vtkOctreePointLocator_h
48
50#include "vtkCommonDataModelModule.h" // For export macro
51
52class vtkCellArray;
53class vtkIdTypeArray;
55class vtkPoints;
56class vtkPolyData;
57
58class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
59{
60public:
62 void PrintSelf(ostream& os, vtkIndent indent) override;
63
65
67
70 vtkSetMacro(MaximumPointsPerRegion, int);
71 vtkGetMacro(MaximumPointsPerRegion, int);
73
75
78 vtkSetMacro(CreateCubicOctants, int);
79 vtkGetMacro(CreateCubicOctants, int);
81
83
89 vtkGetMacro(FudgeFactor, double);
90 vtkSetMacro(FudgeFactor, double);
92
94
98 double* GetBounds() override;
99 void GetBounds(double* bounds) override;
101
103
106 vtkGetMacro(NumberOfLeafNodes, int);
108
112 void GetRegionBounds(int regionID, double bounds[6]);
113
117 void GetRegionDataBounds(int leafNodeID, double bounds[6]);
118
122 int GetRegionContainingPoint(double x, double y, double z);
123
129 void BuildLocator() override;
130
132
136 vtkIdType FindClosestPoint(const double x[3]) override;
137 vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
139
145 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
146
148
153 vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
154 vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
156
161 void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
162
171 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
172
177
181 void FreeSearchStructure() override;
182
187 void GenerateRepresentation(int level, vtkPolyData* pd) override;
188
195 void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
196
197protected:
200
202 vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
203
205
207
211 int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
212 int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
214
216
218
225 vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
226
227 // Recursive helper for public FindPointsWithinRadius
229
230 // Recursive helper for public FindPointsInArea
232
233 // Recursive helper for public FindPointsInArea
235
236 void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
237
238 int DivideTest(int size, int level);
239
241
246 int _FindClosestPointInRegion(int leafNodeId, double x, double y, double z, double& dist2);
247
256 double x, double y, double z, double radius, int skipRegion, double& dist2);
257
259
265
266 double FudgeFactor; // a very small distance, relative to the dataset's size
270
271 float MaxWidth;
272
281
283 void operator=(const vtkOctreePointLocator&) = delete;
284};
285#endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:190
list of point or cell ids
Definition: vtkIdList.h:40
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:43
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.
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
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...
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)
vtkIdType FindClosestPoint(double x, double y, double z, double &dist2)
Return the Id of the point that is closest to the given point.
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:43
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:95
@ level
Definition: vtkX3D.h:401
@ radius
Definition: vtkX3D.h:258
@ size
Definition: vtkX3D.h:259
@ index
Definition: vtkX3D.h:252
int vtkIdType
Definition: vtkType.h:332
#define VTK_NEWINSTANCE