VTK  9.2.6
vtkCellLocator.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkCellLocator.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=========================================================================*/
44
45#ifndef vtkCellLocator_h
46#define vtkCellLocator_h
47
49#include "vtkCommonDataModelModule.h" // For export macro
50#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
51#include "vtkNew.h" // For vtkNew
52
53class vtkIntArray;
54
55class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
56{
57public:
59
63 void PrintSelf(ostream& os, vtkIndent indent) override;
65
71
77
78 // Re-use any superclass signatures that we don't override.
83
90 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
91 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
92
102 int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints* points,
103 vtkIdList* cellIds, vtkGenericCell* cell) override;
104
114 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
115 vtkIdType& cellId, int& subId, double& dist2) override
116 {
117 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
118 }
119
132 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
133 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
134
142 vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* GenCell, int& subId,
143 double pcoords[3], double* weights) override;
144
149 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
150
160 const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
161 {
162 this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
163 }
164
166
169 void FreeSearchStructure() override;
170 void BuildLocator() override;
171 void ForceBuildLocator() override;
172 void GenerateRepresentation(int level, vtkPolyData* pd) override;
174
175 VTK_DEPRECATED_IN_9_2_0("This method is deprecated because LazyEvaluation has been deprecated")
176 virtual void BuildLocatorIfNeeded() {}
177
181 virtual vtkIdList* GetCells(int bucket);
182
187 virtual int GetNumberOfBuckets(void);
188
192 void ShallowCopy(vtkAbstractCellLocator* locator) override;
193
194protected:
196 ~vtkCellLocator() override;
197
198 void BuildLocatorInternal() override;
199
200 //------------------------------------------------------------------------------
202 {
203 public:
205
207
208 inline void Reset();
209
210 inline int* GetPoint(int i);
211
212 inline int InsertNextPoint(int* x);
213
214 protected:
216 };
217
218 void GetOverlappingBuckets(vtkNeighborCells& buckets, const double x[3], double dist,
219 int prevMinLevel[3], int prevMaxLevel[3]);
220
221 inline void GetBucketIndices(const double x[3], int ijk[3]);
222
223 double Distance2ToBucket(const double x[3], int nei[3]);
224 double Distance2ToBounds(const double x[3], double bounds[6]);
225
226 int NumberOfOctants; // number of octants in tree
227 double Bounds[6]; // bounding box root octant
228 double H[3]; // width of leaf octant in x-y-z directions
229 int NumberOfDivisions; // number of "leaf" octant sub-divisions
230 std::shared_ptr<std::vector<vtkSmartPointer<vtkIdList>>> TreeSharedPtr;
232
233 void MarkParents(const vtkSmartPointer<vtkIdList>&, int, int, int, int, int);
234 int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
236 int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
237 void ComputeOctantBounds(double octantBounds[6], int i, int j, int k);
238
239private:
240 vtkCellLocator(const vtkCellLocator&) = delete;
241 void operator=(const vtkCellLocator&) = delete;
242};
243
244#endif
virtual vtkIdType FindCell(double x[3])
Returns the Id of the cell containing the point, returns -1 if no cell found.
virtual void SetNumberOfCellsPerNode(int)
Specify the preferred/maximum number of cells in each node/bucket.
virtual void FindClosestPoint(const double x[3], double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
Return the closest point and the cell which is closest to the point x.
virtual void FindCellsAlongLine(const double p1[3], const double p2[3], double tolerance, vtkIdList *cells)
Take the passed line segment and intersect it with the data set.
virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
Return the closest point within a specified radius and the cell which is closest to the point x.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)
Return intersection point (if any) of finite line with cells contained in cell locator.
object to represent cell connectivity
virtual vtkIdList * GetCells(int bucket)
Get the cells in a particular bucket.
void MarkParents(const vtkSmartPointer< vtkIdList > &, int, int, int, int, int)
void FreeSearchStructure() override
Satisfy vtkLocator abstract interface.
static vtkCellLocator * New()
Construct with automatic computation of divisions, averaging 25 cells per bucket.
~vtkCellLocator() override
double Distance2ToBounds(const double x[3], double bounds[6])
int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType &idx)
int GetNumberOfCellsPerBucket()
void GetBucketIndices(const double x[3], int ijk[3])
void GenerateRepresentation(int level, vtkPolyData *pd) override
Satisfy vtkLocator abstract interface.
void GetOverlappingBuckets(vtkNeighborCells &buckets, const double x[3], double dist, int prevMinLevel[3], int prevMaxLevel[3])
virtual int GetNumberOfBuckets(void)
Return number of buckets available.
void SetNumberOfCellsPerBucket(int N)
Specify the average number of cells in each octant.
void FindCellsWithinBounds(double *bbox, vtkIdList *cells) override
Return a list of unique cell ids inside of a given bounding box.
int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell) override
Take the passed line segment and intersect it with the data set.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to print and obtain type-related information.
void ComputeOctantBounds(double octantBounds[6], int i, int j, int k)
double Distance2ToBucket(const double x[3], int nei[3])
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId, vtkIdType &cellId, vtkGenericCell *cell) override
Return intersection point (if any) AND the cell which was intersected by the finite line.
std::shared_ptr< std::vector< vtkSmartPointer< vtkIdList > > > TreeSharedPtr
void GenerateFace(int face, int numDivs, int i, int j, int k, vtkPoints *pts, vtkCellArray *polys)
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
void FindCellsAlongLine(const double p1[3], const double p2[3], double tolerance, vtkIdList *cellsIds) override
Take the passed line segment and intersect it with the data set.
void ForceBuildLocator() override
Satisfy vtkLocator abstract interface.
void ShallowCopy(vtkAbstractCellLocator *locator) override
Shallow copy of a vtkCellLocator.
vtkSmartPointer< vtkIdList > * Tree
vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2, int &inside) override
Return the closest point within a specified radius and the cell which is closest to the point x.
void BuildLocator() override
Satisfy vtkLocator abstract interface.
void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2) override
Return the closest point and the cell which is closest to the point x.
virtual void BuildLocatorIfNeeded()
vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell *GenCell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:31
a simple class to control print indentation
Definition vtkIndent.h:34
dynamic, self-adjusting array of int
Definition vtkIntArray.h:40
Allocate and hold a VTK object.
Definition vtkNew.h:56
represent and manipulate 3D points
Definition vtkPoints.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:85
Hold a reference to a vtkObjectBase instance.
#define VTK_DEPRECATED_IN_9_2_0(reason)
int vtkIdType
Definition vtkType.h:332