VTK  9.4.20241218
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);
59 vtkGetObjectMacro(CellLocator, vtkAbstractCellLocator);
61
63
70 vtkSetMacro(UseAllPoints, bool);
71 vtkBooleanMacro(UseAllPoints, bool);
72 vtkGetMacro(UseAllPoints, bool);
74
79 vtkGetMacro(IsLinearTransformation, bool);
80
81 // Re-use any superclass signatures that we don't override.
86
93 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
94 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
95
105 int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
106 vtkIdList* cellIds, vtkGenericCell* cell) override;
107
115 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
116 vtkIdType& cellId, int& subId, double& dist2) override
117 {
118 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
119 }
120
131 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
132 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
133
140 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
141
149 const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
150 {
151 this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
152 }
153
161 const double o[3], const double n[3], double tolerance, vtkIdList* cells) override;
162
170 vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* cell, int& subId,
171 double pcoords[3], double* weights) override;
172
177 bool InsideCellBounds(double x[3], vtkIdType cellId) override;
178
180
183 void GenerateRepresentation(int level, vtkPolyData* pd) override;
184 void FreeSearchStructure() override;
185 void BuildLocator() override;
186 void ForceBuildLocator() override;
188
194 void ShallowCopy(vtkAbstractCellLocator* locator) override;
195
196protected:
199
200 void BuildLocatorInternal() override;
201
203
206 bool IsLinearTransformation = false;
207 bool UseAllPoints = false;
208
210
211private:
213 void operator=(const vtkLinearTransformCellLocator&) = delete;
214};
215
216VTK_ABI_NAMESPACE_END
217#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
Cell locator adaptor to perform cell Location on datasets that are a linear transformation of the ori...
void ForceBuildLocator() override
Satisfy vtkLocator abstract interface.
bool InsideCellBounds(double x[3], vtkIdType cellId) 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 tolerance, vtkIdList *cells) override
Take the passed line segment and intersect it with the data set.
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.
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.
vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell *cell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
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 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.
represent and manipulate 3D points
Definition vtkPoints.h:139
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:315