VTK  9.3.20240423
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
40#ifndef vtkStaticCellLocator_h
41#define vtkStaticCellLocator_h
42
44#include "vtkCommonDataModelModule.h" // For export macro
45
46// Forward declarations for PIMPL
47VTK_ABI_NAMESPACE_BEGIN
48struct vtkCellBinner;
49struct vtkCellProcessor;
50
51class VTKCOMMONDATAMODEL_EXPORT vtkStaticCellLocator : public vtkAbstractCellLocator
52{
53 friend struct vtkCellBinner;
54 friend struct vtkCellProcessor;
55
56public:
58
63 void PrintSelf(ostream& os, vtkIndent indent) override;
65
67
73 vtkSetVector3Macro(Divisions, int);
74 vtkGetVectorMacro(Divisions, int, 3);
76
78
92 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
93 vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
95
103 bool GetLargeIds() { return this->LargeIds; }
104
105 // Re-use any superclass signatures that we don't override.
110
117 int IntersectWithLine(const double a0[3], const double a1[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
157 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
158 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
159
164 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
165
175 const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
176 {
177 this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
178 }
179
187 const double o[3], const double n[3], double tolerance, vtkIdList* cells) override;
188
196 vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* GenCell, int& subId,
197 double pcoords[3], double* weights) override;
198
203 bool InsideCellBounds(double x[3], vtkIdType cellId) override;
204
206
209 void GenerateRepresentation(int level, vtkPolyData* pd) override;
210 void FreeSearchStructure() override;
211 void BuildLocator() override;
212 void ForceBuildLocator() override;
214
220 void ShallowCopy(vtkAbstractCellLocator* locator) override;
221
222protected:
225
226 void BuildLocatorInternal() override;
227
228 double Bounds[6]; // Bounding box of the whole dataset
229 int Divisions[3]; // Number of sub-divisions in x-y-z directions
230 double H[3]; // Width of each bin in x-y-z directions
231
232 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
233 bool LargeIds; // indicate whether integer ids are small or large
234
235 // Support PIMPLd implementation
236 vtkCellBinner* Binner; // Does the binning
237 vtkCellProcessor* Processor; // Invokes methods (templated subclasses)
238
239private:
241 void operator=(const vtkStaticCellLocator&) = delete;
242};
243
244VTK_ABI_NAMESPACE_END
245#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 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.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
a simple class to control print indentation
Definition vtkIndent.h:108
represent and manipulate 3D points
Definition vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
perform fast cell location operations
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 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 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.
void FindCellsWithinBounds(double *bbox, vtkIdList *cells) override
Return a list of unique cell ids inside of a given bounding box.
void FindCellsAlongPlane(const double o[3], const double n[3], double tolerance, vtkIdList *cells) override
Take the passed line segment and intersect it with the data set.
void ShallowCopy(vtkAbstractCellLocator *locator) override
Shallow copy of a vtkStaticCellLocator.
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.
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.
bool InsideCellBounds(double x[3], vtkIdType cellId) override
Quickly test if a point is inside the bounds of a particular cell.
~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 IntersectWithLine(const double a0[3], const double a1[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.
int vtkIdType
Definition vtkType.h:315
#define VTK_ID_MAX
Definition vtkType.h:319