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 
129 #ifndef vtkOctreePointLocator_h
130 #define vtkOctreePointLocator_h
131 
132 #include "vtkAbstractPointLocator.h"
133 #include "vtkCommonDataModelModule.h" // For export macro
134 
135 class vtkCellArray;
136 class vtkIdTypeArray;
138 class vtkPoints;
139 class vtkPolyData;
140 
141 class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
142 {
143 public:
145  void PrintSelf(ostream& os, vtkIndent indent) override;
146 
148 
150 
153  vtkSetMacro(MaximumPointsPerRegion, int);
154  vtkGetMacro(MaximumPointsPerRegion, int);
156 
158 
161  vtkSetMacro(CreateCubicOctants, int);
162  vtkGetMacro(CreateCubicOctants, int);
164 
166 
172  vtkGetMacro(FudgeFactor, double);
173  vtkSetMacro(FudgeFactor, double);
175 
177 
181  double* GetBounds() override;
182  void GetBounds(double* bounds) override;
184 
186 
189  vtkGetMacro(NumberOfLeafNodes, int);
191 
195  void GetRegionBounds(int regionID, double bounds[6]);
196 
200  void GetRegionDataBounds(int leafNodeID, double bounds[6]);
201 
205  int GetRegionContainingPoint(double x, double y, double z);
206 
214  void BuildLocator() override;
215 
219  void ForceBuildLocator() override;
220 
222 
226  vtkIdType FindClosestPoint(const double x[3]) override;
227  vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
229 
235  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
236 
238 
243  vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
244  vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
246 
251  void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
252 
261  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
262 
267 
271  void FreeSearchStructure() override;
272 
277  void GenerateRepresentation(int level, vtkPolyData* pd) override;
278 
285  void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
286 
287 protected:
290 
291  void BuildLocatorInternal() override;
292 
294  vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
295 
297 
299 
303  int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
304  int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
306 
308 
310 
317  vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
318 
319  // Recursive helper for public FindPointsWithinRadius
321 
322  // Recursive helper for public FindPointsInArea
324 
325  // Recursive helper for public FindPointsInArea
327 
328  void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
329 
330  int DivideTest(int size, int level);
331 
333 
338  int FindClosestPointInRegion_(int leafNodeId, double x, double y, double z, double& dist2);
339 
348  double x, double y, double z, double radius, int skipRegion, double& dist2);
349 
351 
357 
358  double FudgeFactor; // a very small distance, relative to the dataset's size
362 
363  float MaxWidth;
364 
373 
375  void operator=(const vtkOctreePointLocator&) = delete;
376 };
377 #endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:296
list of point or cell ids
Definition: vtkIdList.h:143
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:119
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:149
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:200
@ 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