VTK  9.2.6
vtkLinearTransformCellLocator.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkLinearTransformCellLocator.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
38#ifndef vtkLinearTransformCellLocator_h
39#define vtkLinearTransformCellLocator_h
40
42#include "vtkFiltersFlowPathsModule.h" // For export macro
43#include "vtkSmartPointer.h" // For vtkSmartPointer
44
45class vtkTransform;
46
47class VTKFILTERSFLOWPATHS_EXPORT vtkLinearTransformCellLocator : public vtkAbstractCellLocator
48{
49public:
51
56 void PrintSelf(ostream& os, vtkIndent indent) override;
58
60
65 virtual void SetCellLocator(vtkAbstractCellLocator* locator);
68
70
77 vtkSetMacro(UseAllPoints, bool);
78 vtkBooleanMacro(UseAllPoints, bool);
79 vtkGetMacro(UseAllPoints, bool);
81
86 vtkGetMacro(IsLinearTransformation, bool);
87
88 // Re-use any superclass signatures that we don't override.
93
100 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
101 double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
102
112 int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints* points,
113 vtkIdList* cellIds, vtkGenericCell* cell) override;
114
122 void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
123 vtkIdType& cellId, int& subId, double& dist2) override
124 {
125 this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
126 }
127
138 vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
139 vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
140
147 void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
148
156 const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
157 {
158 this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
159 }
160
168 const double o[3], const double n[3], double tolerance, vtkIdList* cells) override;
169
177 vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* cell, int& subId,
178 double pcoords[3], double* weights) override;
179
184 bool InsideCellBounds(double x[3], vtkIdType cellId) override;
185
187
190 void GenerateRepresentation(int level, vtkPolyData* pd) override;
191 void FreeSearchStructure() override;
192 void BuildLocator() override;
193 void ForceBuildLocator() override;
195
199 void ShallowCopy(vtkAbstractCellLocator* locator) override;
200
201protected:
204
205 void BuildLocatorInternal() override;
206
208
212 bool UseAllPoints = false;
213
215
216private:
218 void operator=(const vtkLinearTransformCellLocator&) = delete;
219};
220
221#endif
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 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 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:31
a simple class to control print indentation
Definition vtkIndent.h:34
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
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.
int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell) override
Take the passed line segment and intersect it with the data set.
represent and manipulate 3D points
Definition vtkPoints.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:85
Hold a reference to a vtkObjectBase instance.
describes linear transformations via a 4x4 matrix
int vtkIdType
Definition vtkType.h:332