VTK  9.6.20260518
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
113
114#ifndef vtkCellLocator_h
115#define vtkCellLocator_h
116
118#include "vtkCommonDataModelModule.h" // For export macro
119#include "vtkNew.h" // For vtkNew
120
121VTK_ABI_NAMESPACE_BEGIN
122class vtkIntArray;
123
124class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
125{
126public:
128
132 void PrintSelf(ostream& os, vtkIndent indent) override;
134
140
146
147 // Reuse any superclass signatures that we don't override.
153
160 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
161 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
162
172 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
173 vtkIdList* cellIds, vtkGenericCell* cell) override;
174
184 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
185 vtkIdType& cellId, int& subId, double& dist2) override
186 {
187 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
188 }
189
210 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
211 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
212
220 vtkIdType FindCell(double x[3], double tol2, vtkGenericCell* GenCell, int& subId,
221 double pcoords[3], double* weights) override;
222
227 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
228
230
233 void FreeSearchStructure() override;
234 void BuildLocator() override;
235 void ForceBuildLocator() override;
236 void GenerateRepresentation(int level, vtkPolyData* pd) override;
238
242 virtual vtkIdList* GetCells(int bucket);
243
248 virtual int GetNumberOfBuckets();
249
255 void ShallowCopy(vtkAbstractCellLocator* locator) override;
256
257protected:
259 ~vtkCellLocator() override;
260
261 void BuildLocatorInternal() override;
262
263 //------------------------------------------------------------------------------
265 {
266 public:
268
270
271 inline void Reset();
272
273 inline int* GetPoint(int i);
274
275 inline int InsertNextPoint(int* x);
276
277 protected:
279 };
280
281 void GetOverlappingBuckets(vtkNeighborCells& buckets, const double x[3], double dist,
282 int prevMinLevel[3], int prevMaxLevel[3]);
283
284 inline void GetBucketIndices(const double x[3], int ijk[3]);
285
286 double Distance2ToBucket(const double x[3], int nei[3]);
287 double Distance2ToBounds(const double x[3], double bounds[6]);
288
289 int NumberOfOctants; // number of octants in tree
290 double Bounds[6]; // bounding box root octant
291 double H[3]; // width of leaf octant in x-y-z directions
292 int NumberOfDivisions; // number of "leaf" octant sub-divisions
293 std::shared_ptr<std::vector<vtkSmartPointer<vtkIdList>>> TreeSharedPtr;
295
296 void MarkParents(const vtkSmartPointer<vtkIdList>&, int, int, int, int, int);
297 int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
299 int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
300 void ComputeOctantBounds(double octantBounds[6], int i, int j, int k);
301
302private:
303 vtkCellLocator(const vtkCellLocator&) = delete;
304 void operator=(const vtkCellLocator&) = delete;
305};
306
307VTK_ABI_NAMESPACE_END
308#endif
virtual void SetNumberOfCellsPerNode(int)
Specify the preferred/maximum number of cells in each node/bucket.
vtkIdType FindCell(double x[3])
Returns the Id of the cell containing the point, returns -1 if no cell found.
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.
void FindCellsAlongLine(const double p1[3], const double p2[3], double tol, vtkIdList *cells)
Take the passed line segment and intersect it with the data set.
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.
vtkAbstractCellLocator()
Find the cell containing a given point.
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])
vtkIdType FindCell(double x[3], double tol2, vtkGenericCell *GenCell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
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 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.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:135
a simple class to control print indentation
Definition vtkIndent.h:108
dynamic, self-adjusting array of int
Allocate and hold a VTK object.
Definition vtkNew.h:168
represent and manipulate 3D points
Definition vtkPoints.h:140
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
int vtkIdType
Definition vtkType.h:363