VTK  9.3.0
vtkCellLocator.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
36#ifndef vtkCellLocator_h
37#define vtkCellLocator_h
38
40#include "vtkCommonDataModelModule.h" // For export macro
41#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
42#include "vtkNew.h" // For vtkNew
43
44VTK_ABI_NAMESPACE_BEGIN
45class vtkIntArray;
46
47class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
48{
49public:
51
55 void PrintSelf(ostream& os, vtkIndent indent) override;
57
63
69
70 // Re-use any superclass signatures that we don't override.
75
82 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
83 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
84
94 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
95 vtkIdList* cellIds, vtkGenericCell* cell) override;
96
106 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
107 vtkIdType& cellId, int& subId, double& dist2) override
108 {
109 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
110 }
111
124 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
125 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
126
134 vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* GenCell, int& subId,
135 double pcoords[3], double* weights) override;
136
141 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
142
152 const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
153 {
154 this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
155 }
156
158
161 void FreeSearchStructure() override;
162 void BuildLocator() override;
163 void ForceBuildLocator() override;
164 void GenerateRepresentation(int level, vtkPolyData* pd) override;
166
167 VTK_DEPRECATED_IN_9_2_0("This method is deprecated because LazyEvaluation has been deprecated")
168 virtual void BuildLocatorIfNeeded() {}
169
173 virtual vtkIdList* GetCells(int bucket);
174
179 virtual int GetNumberOfBuckets();
180
184 void ShallowCopy(vtkAbstractCellLocator* locator) override;
185
186protected:
188 ~vtkCellLocator() override;
189
190 void BuildLocatorInternal() override;
191
192 //------------------------------------------------------------------------------
194 {
195 public:
197
199
200 inline void Reset();
201
202 inline int* GetPoint(int i);
203
204 inline int InsertNextPoint(int* x);
205
206 protected:
208 };
209
210 void GetOverlappingBuckets(vtkNeighborCells& buckets, const double x[3], double dist,
211 int prevMinLevel[3], int prevMaxLevel[3]);
212
213 inline void GetBucketIndices(const double x[3], int ijk[3]);
214
215 double Distance2ToBucket(const double x[3], int nei[3]);
216 double Distance2ToBounds(const double x[3], double bounds[6]);
217
218 int NumberOfOctants; // number of octants in tree
219 double Bounds[6]; // bounding box root octant
220 double H[3]; // width of leaf octant in x-y-z directions
221 int NumberOfDivisions; // number of "leaf" octant sub-divisions
222 std::shared_ptr<std::vector<vtkSmartPointer<vtkIdList>>> TreeSharedPtr;
224
225 void MarkParents(const vtkSmartPointer<vtkIdList>&, int, int, int, int, int);
226 int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
228 int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
229 void ComputeOctantBounds(double octantBounds[6], int i, int j, int k);
230
231private:
232 vtkCellLocator(const vtkCellLocator&) = delete;
233 void operator=(const vtkCellLocator&) = delete;
234};
235
236VTK_ABI_NAMESPACE_END
237#endif
an abstract base class for locators which find cells
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 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
octree-based spatial search object to quickly locate cells
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()
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.
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)
int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell) override
Take the passed line segment and intersect it with the data set.
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.
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:23
a simple class to control print indentation
Definition vtkIndent.h:29
dynamic, self-adjusting array of int
Definition vtkIntArray.h:35
Allocate and hold a VTK object.
Definition vtkNew.h:51
represent and manipulate 3D points
Definition vtkPoints.h:29
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
Hold a reference to a vtkObjectBase instance.
#define VTK_DEPRECATED_IN_9_2_0(reason)
int vtkIdType
Definition vtkType.h:315