VTK  9.6.20260516
vtkStaticCellLocator.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
38
39#ifndef vtkStaticCellLocator_h
40#define vtkStaticCellLocator_h
41
43#include "vtkCommonDataModelModule.h" // For export macro
44
45// Forward declarations for PIMPL
46VTK_ABI_NAMESPACE_BEGIN
47struct vtkCellBinner;
48struct vtkCellProcessor;
49
50class VTKCOMMONDATAMODEL_EXPORT vtkStaticCellLocator : public vtkAbstractCellLocator
51{
52 friend struct vtkCellBinner;
53 friend struct vtkCellProcessor;
54
55public:
57
62 void PrintSelf(ostream& os, vtkIndent indent) override;
64
66
72 vtkSetVector3Macro(Divisions, int);
73 vtkGetVectorMacro(Divisions, int, 3);
75
77
91 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
94
102 bool GetLargeIds() { return this->LargeIds; }
103
104 // Reuse any superclass signatures that we don't override.
110
117 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
118 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
119
129 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
130 vtkIdList* cellIds, vtkGenericCell* cell) override;
131
141 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
142 vtkIdType& cellId, int& subId, double& dist2) override
143 {
144 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
145 }
146
167 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
168 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
169
174 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
175
183 const double o[3], const double n[3], double tol, vtkIdList* cells) override;
184
192 vtkIdType FindCell(double x[3], double tol2, vtkGenericCell* GenCell, int& subId,
193 double pcoords[3], double* weights) override;
194
199 bool InsideCellBounds(double x[3], vtkIdType cellId, double tol = 0.0) override;
200
202
205 void GenerateRepresentation(int level, vtkPolyData* pd) override;
206 void FreeSearchStructure() override;
207 void BuildLocator() override;
208 void ForceBuildLocator() override;
210
216 void ShallowCopy(vtkAbstractCellLocator* locator) override;
217
218protected:
221
222 void BuildLocatorInternal() override;
223
224 double Bounds[6]; // Bounding box of the whole dataset
225 int Divisions[3]; // Number of sub-divisions in x-y-z directions
226 double H[3]; // Width of each bin in x-y-z directions
227
228 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
229 bool LargeIds; // indicate whether integer ids are small or large
230
231 // Support PIMPLd implementation
232 vtkCellBinner* Binner; // Does the binning
233 vtkCellProcessor* Processor; // Invokes methods (templated subclasses)
234
235private:
237 void operator=(const vtkStaticCellLocator&) = delete;
238};
239
240VTK_ABI_NAMESPACE_END
241#endif
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.
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
represent and manipulate 3D points
Definition vtkPoints.h:140
concrete dataset represents vertices, lines, polygons, and triangle strips
static vtkStaticCellLocator * New()
Standard methods to instantiate, print and obtain type-related information.
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 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 tol2, vtkGenericCell *GenCell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
void FindCellsWithinBounds(double *bbox, vtkIdList *cells) override
Return a list of unique cell ids inside of a given bounding box.
void ShallowCopy(vtkAbstractCellLocator *locator) override
Shallow copy of a vtkStaticCellLocator.
bool InsideCellBounds(double x[3], vtkIdType cellId, double tol=0.0) override
Quickly test if a point is inside the bounds of a particular cell.
void FindCellsAlongPlane(const double o[3], const double n[3], double tol, vtkIdList *cells) override
Take the passed line segment and intersect it with the data set.
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 BuildLocator() override
Satisfy vtkLocator abstract interface.
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.
~vtkStaticCellLocator() override
vtkCellProcessor * Processor
bool GetLargeIds()
Inform the user as to whether large ids are being used.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to instantiate, print and obtain type-related information.
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
void GenerateRepresentation(int level, vtkPolyData *pd) override
Satisfy vtkLocator abstract interface.
void FreeSearchStructure() override
Satisfy vtkLocator abstract interface.
void ForceBuildLocator() override
Satisfy vtkLocator abstract interface.
int vtkIdType
Definition vtkType.h:363
#define VTK_ID_MAX
Definition vtkType.h:367