VTK  9.4.20241108
vtkAbstractCellLocator.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
28#ifndef vtkAbstractCellLocator_h
29#define vtkAbstractCellLocator_h
30
31#include "vtkCommonDataModelModule.h" // For export macro
32#include "vtkLocator.h"
33#include "vtkNew.h" // For vtkNew
34
35#include <memory> // For shared_ptr
36#include <vector> // For Weights
37
38VTK_ABI_NAMESPACE_BEGIN
39class vtkCellArray;
40class vtkGenericCell;
41class vtkIdList;
42class vtkPoints;
43
44class VTKCOMMONDATAMODEL_EXPORT vtkAbstractCellLocator : public vtkLocator
45{
46public:
48 void PrintSelf(ostream& os, vtkIndent indent) override;
49
51
57 vtkSetClampMacro(NumberOfCellsPerNode, int, 1, VTK_INT_MAX);
58 vtkGetMacro(NumberOfCellsPerNode, int);
60
62
69 vtkSetMacro(CacheCellBounds, vtkTypeBool);
70 vtkGetMacro(CacheCellBounds, vtkTypeBool);
71 vtkBooleanMacro(CacheCellBounds, vtkTypeBool);
73
79
81
86 vtkSetMacro(RetainCellLists, vtkTypeBool);
87 vtkGetMacro(RetainCellLists, vtkTypeBool);
88 vtkBooleanMacro(RetainCellLists, vtkTypeBool);
90
97 virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
98 double x[3], double pcoords[3], int& subId);
99
106 virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
107 double x[3], double pcoords[3], int& subId, vtkIdType& cellId);
108
117 virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
118 double x[3], double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell);
119
135 virtual int IntersectWithLine(
136 const double p1[3], const double p2[3], vtkPoints* points, vtkIdList* cellIds);
137
147 virtual int IntersectWithLine(
148 const double p1[3], const double p2[3], double tol, vtkPoints* points, vtkIdList* cellIds);
149
161 virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol,
162 vtkPoints* points, vtkIdList* cellIds, vtkGenericCell* cell);
163
174 virtual void FindClosestPoint(
175 const double x[3], double closestPoint[3], vtkIdType& cellId, int& subId, double& dist2);
176
189 virtual void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
190 vtkIdType& cellId, int& subId, double& dist2);
191
202 virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
203 vtkIdType& cellId, int& subId, double& dist2);
204
217 virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
218 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2);
219
234 virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
235 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside);
236
243 virtual void FindCellsWithinBounds(double* bbox, vtkIdList* cells);
244
256 virtual void FindCellsAlongLine(
257 const double p1[3], const double p2[3], double tolerance, vtkIdList* cells);
258
269 const double o[3], const double n[3], double tolerance, vtkIdList* cells);
270
277 virtual vtkIdType FindCell(double x[3]);
278
280
288 double x[3], double tol2, vtkGenericCell* GenCell, double pcoords[3], double* weights);
289 virtual vtkIdType FindCell(double x[3], double tol2, vtkGenericCell* GenCell, int& subId,
290 double pcoords[3], double* weights);
292
298 virtual bool InsideCellBounds(double x[3], vtkIdType cell_ID);
299
306
307protected:
310
312
319 virtual bool StoreCellBounds();
320 virtual void FreeCellBounds();
322
328
333 std::shared_ptr<std::vector<double>> CellBoundsSharedPtr;
334 double* CellBounds; // The is just used for simplicity in the internal code
335
340
341 static bool IsInBounds(const double bounds[6], const double x[3], double tol = 0.0);
342
343 /*
344 * This function should be used ONLY after the locator is built.
345 * cellBoundsPtr should be assigned to a double cellBounds[6] BEFORE calling this function.
346 */
347 void GetCellBounds(vtkIdType cellId, double*& cellBoundsPtr);
348
355 std::vector<double> Weights;
356
357private:
359 void operator=(const vtkAbstractCellLocator&) = delete;
360};
361
362VTK_ABI_NAMESPACE_END
363#endif
an abstract base class for locators which find cells
virtual int IntersectWithLine(const double p1[3], const double p2[3], vtkPoints *points, vtkIdList *cellIds)
Take the passed line segment and intersect it with the data set.
virtual vtkIdType FindCell(double x[3])
Returns the Id of the cell containing the point, returns -1 if no cell found.
void UpdateInternalWeights()
To be called in FindCell(double[3]).
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void FindCellsAlongPlane(const double o[3], const double n[3], double tolerance, vtkIdList *cells)
Given an unbounded plane defined by an origin o[3] and unit normal n[3], return the list of unique ce...
vtkTimeStamp WeightsTime
This time stamp helps us decide if we want to update internal Weights array size.
virtual void ShallowCopy(vtkAbstractCellLocator *)
Shallow copy of a vtkAbstractCellLocator.
~vtkAbstractCellLocator() override
vtkNew< vtkGenericCell > GenericCell
virtual vtkIdType FindCell(double x[3], double tol2, vtkGenericCell *GenCell, double pcoords[3], double *weights)
Find the cell containing a given point.
virtual bool StoreCellBounds()
This command is used internally by the locator to copy all cell Bounds into the internal CellBounds a...
virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, 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 GetCellBounds(vtkIdType cellId, double *&cellBoundsPtr)
virtual void FindCellsWithinBounds(double *bbox, vtkIdList *cells)
Return a list of unique cell ids inside of a given bounding box.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints *points, vtkIdList *cellIds)
Take the passed line segment and intersect it with the data set.
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 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell)
Take the passed line segment and intersect it with the data set.
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 void FreeCellBounds()
This command is used internally by the locator to copy all cell Bounds into the internal CellBounds a...
virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2, int &inside)
Return the closest point within a specified radius and the cell which is closest to the point x.
std::vector< double > Weights
This array is resized so that it can fit points from the cell hosting the most in the input data set.
static bool IsInBounds(const double bounds[6], const double x[3], double tol=0.0)
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 ComputeCellBounds()
This function can be used either internally or externally to compute only the cached cell bounds if C...
std::shared_ptr< std::vector< double > > CellBoundsSharedPtr
virtual bool InsideCellBounds(double x[3], vtkIdType cell_ID)
Quickly test if a point is inside the bounds of a particular cell.
virtual void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2)
Return the closest point and the cell which is closest to the point x.
virtual vtkIdType FindCell(double x[3], double tol2, vtkGenericCell *GenCell, int &subId, double pcoords[3], double *weights)
Find the cell containing a given point.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId, vtkIdType &cellId)
Return intersection point (if any) AND the cell which was intersected by the finite line.
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.
virtual 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)
Return intersection point (if any) AND the cell which was intersected by the finite line.
object to represent cell connectivity
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
abstract base class for objects that accelerate spatial searches
Definition vtkLocator.h:78
Allocate and hold a VTK object.
Definition vtkNew.h:167
represent and manipulate 3D points
Definition vtkPoints.h:139
record modification and/or execution time
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
#define VTK_INT_MAX
Definition vtkType.h:144