VTK  9.6.20260516
vtkLinearTransformCellLocator.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
29#ifndef vtkLinearTransformCellLocator_h
30#define vtkLinearTransformCellLocator_h
31
33#include "vtkFiltersFlowPathsModule.h" // For export macro
34#include "vtkSmartPointer.h" // For vtkSmartPointer
35
36VTK_ABI_NAMESPACE_BEGIN
37
38class vtkTransform;
39
40class VTKFILTERSFLOWPATHS_EXPORT vtkLinearTransformCellLocator : public vtkAbstractCellLocator
41{
42public:
44
49 void PrintSelf(ostream& os, vtkIndent indent) override;
51
53
58 virtual void SetCellLocator(vtkAbstractCellLocator* locator);
61
63
70 vtkSetMacro(UseAllPoints, bool);
71 vtkBooleanMacro(UseAllPoints, bool);
72 vtkGetMacro(UseAllPoints, bool);
74
79 vtkGetMacro(IsLinearTransformation, bool);
80
81 // Reuse any superclass signatures that we don't override.
87
94 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
95 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
96
106 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
107 vtkIdList* cellIds, vtkGenericCell* cell) override;
108
116 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
117 vtkIdType& cellId, int& subId, double& dist2) override
118 {
119 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
120 }
121
142 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
143 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
144
151 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
152
160 const double o[3], const double n[3], double tol, vtkIdList* cells) override;
161
169 vtkIdType FindCell(double x[3], double tol2, vtkGenericCell* cell, int& subId, double pcoords[3],
170 double* weights) override;
171
176 bool InsideCellBounds(double x[3], vtkIdType cellId, double tol = 0.0) override;
177
179
182 void GenerateRepresentation(int level, vtkPolyData* pd) override;
183 void FreeSearchStructure() override;
184 void BuildLocator() override;
185 void ForceBuildLocator() override;
187
193 void ShallowCopy(vtkAbstractCellLocator* locator) override;
194
195protected:
198
199 void BuildLocatorInternal() override;
200
202
206 bool UseAllPoints = false;
207
209
210private:
212 void operator=(const vtkLinearTransformCellLocator&) = delete;
213};
214
215VTK_ABI_NAMESPACE_END
216#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
void ForceBuildLocator() override
Satisfy vtkLocator abstract interface.
vtkIdType FindCell(double x[3], double tol2, vtkGenericCell *cell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
vtkSmartPointer< vtkTransform > Transform
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to instantiate, print and obtain type-related information.
void GenerateRepresentation(int level, vtkPolyData *pd) override
Satisfy vtkLocator abstract interface.
void FindCellsWithinBounds(double *bbox, vtkIdList *cells) override
Return a list of unique cell ids inside of a given bounding box.
static vtkLinearTransformCellLocator * New()
Standard methods to instantiate, print and obtain type-related information.
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.
vtkSmartPointer< vtkTransform > InverseTransform
~vtkLinearTransformCellLocator() override
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 ShallowCopy(vtkAbstractCellLocator *locator) override
Shallow copy of a vtkLinearTransformCellLocator.
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 FreeSearchStructure() override
Satisfy vtkLocator abstract interface.
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.
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
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.
virtual void SetCellLocator(vtkAbstractCellLocator *locator)
Set/Get the cell locator to be used internally.
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.
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.
describes linear transformations via a 4x4 matrix
int vtkIdType
Definition vtkType.h:363