VTK  9.0.1
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 
37 #ifndef vtkOctreePointLocator_h
38 #define vtkOctreePointLocator_h
39 
41 #include "vtkCommonDataModelModule.h" // For export macro
42 
43 class vtkCellArray;
44 class vtkIdTypeArray;
46 class vtkPoints;
47 class vtkPolyData;
48 
49 class VTKCOMMONDATAMODEL_EXPORT vtkOctreePointLocator : public vtkAbstractPointLocator
50 {
51 public:
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
56 
58 
61  vtkSetMacro(MaximumPointsPerRegion, int);
62  vtkGetMacro(MaximumPointsPerRegion, int);
64 
66 
69  vtkSetMacro(CreateCubicOctants, int);
70  vtkGetMacro(CreateCubicOctants, int);
72 
74 
80  vtkGetMacro(FudgeFactor, double);
81  vtkSetMacro(FudgeFactor, double);
83 
85 
89  double* GetBounds() override;
90  void GetBounds(double* bounds) override;
92 
94 
97  vtkGetMacro(NumberOfLeafNodes, int);
99 
103  void GetRegionBounds(int regionID, double bounds[6]);
104 
108  void GetRegionDataBounds(int leafNodeID, double bounds[6]);
109 
113  int GetRegionContainingPoint(double x, double y, double z);
114 
120  void BuildLocator() override;
121 
123 
127  vtkIdType FindClosestPoint(const double x[3]) override;
128  vtkIdType FindClosestPoint(double x, double y, double z, double& dist2);
130 
136  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
137 
139 
144  vtkIdType FindClosestPointInRegion(int regionId, double* x, double& dist2);
145  vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double& dist2);
147 
152  void FindPointsWithinRadius(double radius, const double x[3], vtkIdList* result) override;
153 
162  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
163 
168 
172  void FreeSearchStructure() override;
173 
178  void GenerateRepresentation(int level, vtkPolyData* pd) override;
179 
186  void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
187 
188 protected:
191 
193  vtkOctreePointLocatorNode** LeafNodeList; // indexed by region/node ID
194 
196 
198 
202  int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
203  int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
205 
207 
209 
216  vtkOctreePointLocatorNode* node, double radiusSquared, const double x[3], vtkIdList* ids);
217 
218  // Recursive helper for public FindPointsWithinRadius
220 
221  // Recursive helper for public FindPointsInArea
223 
224  // Recursive helper for public FindPointsInArea
226 
227  void DivideRegion(vtkOctreePointLocatorNode* node, int* ordering, int level);
228 
229  int DivideTest(int size, int level);
230 
232 
237  int _FindClosestPointInRegion(int leafNodeId, double x, double y, double z, double& dist2);
238 
247  double x, double y, double z, double radius, int skipRegion, double& dist2);
248 
250 
256 
257  double FudgeFactor; // a very small distance, relative to the dataset's size
261 
262  float MaxWidth;
263 
272 
274  void operator=(const vtkOctreePointLocator&) = delete;
275 };
276 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkOctreePointLocator
an octree spatial decomposition of a set of points
Definition: vtkOctreePointLocator.h:50
vtkOctreePointLocator::BuildLocator
void BuildLocator() override
Create the octree decomposition of the cells of the data set or data sets.
vtkOctreePointLocator::FindClosestPointInSphere
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.
vtkOctreePointLocator::FindClosestPoint
vtkIdType FindClosestPoint(double x, double y, double z, double &dist2)
vtkOctreePointLocator::operator=
void operator=(const vtkOctreePointLocator &)=delete
vtkOctreePointLocator::LocatorIds
int * LocatorIds
Definition: vtkOctreePointLocator.h:260
vtkOctreePointLocator::~vtkOctreePointLocator
~vtkOctreePointLocator() override
vtkOctreePointLocator::MaxWidth
float MaxWidth
Definition: vtkOctreePointLocator.h:262
vtkOctreePointLocator::FindPointsInArea
void FindPointsInArea(double *area, vtkIdTypeArray *ids, bool clearArray=true)
Fill ids with points found in area.
vtkOctreePointLocator::GenerateRepresentation
void GenerateRepresentation(int level, vtkPolyData *pd) override
Create a polydata representation of the boundaries of the octree regions.
vtkAbstractPointLocator.h
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkOctreePointLocator::FudgeFactor
double FudgeFactor
Definition: vtkOctreePointLocator.h:257
vtkOctreePointLocator::FreeSearchStructure
void FreeSearchStructure() override
Delete the octree data structure.
vtkOctreePointLocator::AddAllPointsInRegion
void AddAllPointsInRegion(vtkOctreePointLocatorNode *node, vtkIdList *ids)
vtkOctreePointLocator::GetRegionDataBounds
void GetRegionDataBounds(int leafNodeID, double bounds[6])
Get the bounds of the data within the leaf node.
vtkOctreePointLocator::FindClosestNPoints
void FindClosestNPoints(int N, const double x[3], vtkIdList *result) override
Find the closest N points to a position.
vtkOctreePointLocator::_FindClosestPointInRegion
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...
vtkOctreePointLocator::NumberOfLocatorPoints
int NumberOfLocatorPoints
Definition: vtkOctreePointLocator.h:258
vtkOctreePointLocator::FindPointsInArea
void FindPointsInArea(vtkOctreePointLocatorNode *node, double *area, vtkIdTypeArray *ids)
vtkOctreePointLocator::NumberOfLeafNodes
int NumberOfLeafNodes
Definition: vtkOctreePointLocator.h:254
vtkOctreePointLocatorNode
Octree node that has 8 children each of equal size.
Definition: vtkOctreePointLocatorNode.h:47
vtkOctreePointLocator::vtkOctreePointLocator
vtkOctreePointLocator(const vtkOctreePointLocator &)=delete
vtkOctreePointLocator::DivideTest
int DivideTest(int size, int level)
vtkOctreePointLocator::LeafNodeList
vtkOctreePointLocatorNode ** LeafNodeList
Definition: vtkOctreePointLocator.h:193
vtkOctreePointLocator::FindPointsWithinRadius
void FindPointsWithinRadius(vtkOctreePointLocatorNode *node, double radiusSquared, const double x[3], vtkIdList *ids)
Recursive helper for public FindPointsWithinRadius.
vtkOctreePointLocator::New
static vtkOctreePointLocator * New()
vtkOctreePointLocator::Top
vtkOctreePointLocatorNode * Top
Definition: vtkOctreePointLocator.h:192
vtkX3D::level
@ level
Definition: vtkX3D.h:401
vtkOctreePointLocator::CreateCubicOctants
int CreateCubicOctants
If CreateCubicOctants is non-zero, the bounding box of the points will be expanded such that all octa...
Definition: vtkOctreePointLocator.h:271
vtkOctreePointLocator::DivideRegion
void DivideRegion(vtkOctreePointLocatorNode *node, int *ordering, int level)
vtkOctreePointLocator::SetDataBoundsToSpatialBounds
static void SetDataBoundsToSpatialBounds(vtkOctreePointLocatorNode *node)
vtkOctreePointLocator::FindClosestPointInRegion
vtkIdType FindClosestPointInRegion(int regionId, double x, double y, double z, double &dist2)
vtkOctreePointLocator::FindClosestPointWithinRadius
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.
vtkOctreePointLocator::AddAllPointsInRegion
void AddAllPointsInRegion(vtkOctreePointLocatorNode *node, vtkIdTypeArray *ids)
vtkOctreePointLocator::GetRegionBounds
void GetRegionBounds(int regionID, double bounds[6])
Get the spatial bounds of octree region.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:180
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:31
vtkX3D::size
@ size
Definition: vtkX3D.h:259
vtkOctreePointLocator::DeleteAllDescendants
static void DeleteAllDescendants(vtkOctreePointLocatorNode *octant)
vtkOctreePointLocator::GetBounds
void GetBounds(double *bounds) override
vtkOctreePointLocator::AddPolys
void AddPolys(vtkOctreePointLocatorNode *node, vtkPoints *pts, vtkCellArray *polys)
vtkOctreePointLocator::GetRegionContainingPoint
int GetRegionContainingPoint(double x, double y, double z)
Get the id of the leaf region containing the specified location.
vtkOctreePointLocator::GetPointsInRegion
vtkIdTypeArray * GetPointsInRegion(int leafNodeId)
Get a list of the original IDs of all points in a leaf node.
vtkOctreePointLocator::LocatorPoints
float * LocatorPoints
Definition: vtkOctreePointLocator.h:259
vtkIdTypeArray
dynamic, self-adjusting array of vtkIdType
Definition: vtkIdTypeArray.h:36
vtkOctreePointLocator::FindClosestPointInRegion
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.
vtkOctreePointLocator::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkOctreePointLocator::vtkOctreePointLocator
vtkOctreePointLocator()
vtkPolyData
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
vtkOctreePointLocator::FindRegion
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.
vtkAbstractPointLocator
abstract class to quickly locate points in 3-space
Definition: vtkAbstractPointLocator.h:39
vtkOctreePointLocator::BuildLeafNodeList
void BuildLeafNodeList(vtkOctreePointLocatorNode *node, int &index)
vtkOctreePointLocator::FindClosestPoint
vtkIdType FindClosestPoint(const double x[3]) override
Return the Id of the point that is closest to the given point.
vtkX3D::radius
@ radius
Definition: vtkX3D.h:258
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkOctreePointLocator::FindRegion
int FindRegion(vtkOctreePointLocatorNode *node, double x, double y, double z)
vtkOctreePointLocator::FindPointsWithinRadius
void FindPointsWithinRadius(double radius, const double x[3], vtkIdList *result) override
Find all points within a specified radius of position x.
vtkOctreePointLocator::MaximumPointsPerRegion
int MaximumPointsPerRegion
The maximum number of points in a region/octant before it is subdivided.
Definition: vtkOctreePointLocator.h:253
vtkOctreePointLocator::GetBounds
double * GetBounds() override
Get the spatial bounds of the entire octree space.