VTK
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
43#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) VTK_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() VTK_OVERRIDE;
93 void GetBounds(double *bounds) VTK_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
123 void BuildLocator() VTK_OVERRIDE;
124
126
130 vtkIdType FindClosestPoint(const double x[3]) VTK_OVERRIDE;
131 vtkIdType FindClosestPoint(double x, double y, double z, double &dist2);
133
139 vtkIdType FindClosestPointWithinRadius(
140 double radius, const double x[3], double& dist2) VTK_OVERRIDE;
141
143
148 vtkIdType FindClosestPointInRegion(int regionId, double *x, double &dist2);
149 vtkIdType FindClosestPointInRegion(int regionId, double x, double y,
150 double z, double &dist2);
152
157 void FindPointsWithinRadius(
158 double radius, const double x[3], vtkIdList *result) VTK_OVERRIDE;
159
168 void FindClosestNPoints(int N, const double x[3],
169 vtkIdList *result) VTK_OVERRIDE;
170
174 vtkIdTypeArray *GetPointsInRegion(int leafNodeId);
175
179 void FreeSearchStructure() VTK_OVERRIDE;
180
185 void GenerateRepresentation(int level, vtkPolyData *pd) VTK_OVERRIDE;
186
193 void FindPointsInArea(double* area, vtkIdTypeArray* ids, bool clearArray = true);
194
195protected:
196
198 ~vtkOctreePointLocator() VTK_OVERRIDE;
199
201 vtkOctreePointLocatorNode **LeafNodeList; // indexed by region/node ID
202
203 void BuildLeafNodeList(vtkOctreePointLocatorNode* node, int & index);
204
206
210 int FindRegion(vtkOctreePointLocatorNode* node, float x, float y, float z);
211 int FindRegion(vtkOctreePointLocatorNode* node, double x, double y, double z);
213
214 static void SetDataBoundsToSpatialBounds(vtkOctreePointLocatorNode *node);
215
216 static void DeleteAllDescendants(vtkOctreePointLocatorNode* octant);
217
223 void FindPointsWithinRadius(vtkOctreePointLocatorNode* node, double radiusSquared,
224 const double x[3], vtkIdList* ids);
225
226 // Recursive helper for public FindPointsWithinRadius
227 void AddAllPointsInRegion(vtkOctreePointLocatorNode* node, vtkIdList* ids);
228
229 // Recursive helper for public FindPointsInArea
230 void FindPointsInArea(vtkOctreePointLocatorNode* node, double* area, vtkIdTypeArray* ids);
231
232 // Recursive helper for public FindPointsInArea
233 void AddAllPointsInRegion(vtkOctreePointLocatorNode* node, vtkIdTypeArray* ids);
234
235 void DivideRegion(vtkOctreePointLocatorNode *node, int* ordering, int level);
236
237 int DivideTest(int size, int level);
238
239 void AddPolys(vtkOctreePointLocatorNode *node, vtkPoints *pts, vtkCellArray *polys);
240
245 int _FindClosestPointInRegion(int leafNodeId, double x, double y,
246 double z, double &dist2);
247
255 int FindClosestPointInSphere(double x, double y, double z, double radius,
256 int skipRegion, double &dist2);
257
259
262 int MaximumPointsPerRegion;
263 int NumberOfLeafNodes;
265
266 double FudgeFactor; // a very small distance, relative to the dataset's size
267 int NumberOfLocatorPoints;
268 float *LocatorPoints;
269 int *LocatorIds;
270
271 float MaxWidth;
272
280 int CreateCubicOctants;
281
282 vtkOctreePointLocator(const vtkOctreePointLocator&) VTK_DELETE_FUNCTION;
283 void operator=(const vtkOctreePointLocator&) VTK_DELETE_FUNCTION;
284};
285#endif
abstract class to quickly locate points in 3-space
object to represent cell connectivity
Definition: vtkCellArray.h:51
list of point or cell ids
Definition: vtkIdList.h:37
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 PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static vtkOctreePointLocator * New()
double * GetBounds() override
Get the spatial bounds of the entire octree space.
represent and manipulate 3D points
Definition: vtkPoints.h:40
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
@ level
Definition: vtkX3D.h:395
@ radius
Definition: vtkX3D.h:252
@ size
Definition: vtkX3D.h:253
@ index
Definition: vtkX3D.h:246
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
int vtkIdType
Definition: vtkType.h:287