VTK  9.3.0
vtkCubicLine.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) EDF - www.edf.fr
3 // SPDX-License-Identifier: BSD-3-Clause
24 #ifndef vtkCubicLine_h
25 #define vtkCubicLine_h
26 
27 #include "vtkCommonDataModelModule.h" // For export macro
28 #include "vtkNonLinearCell.h"
29 
30 VTK_ABI_NAMESPACE_BEGIN
31 class vtkLine;
32 class vtkDoubleArray;
33 
34 class VTKCOMMONDATAMODEL_EXPORT vtkCubicLine : public vtkNonLinearCell
35 {
36 public:
37  static vtkCubicLine* New();
39  void PrintSelf(ostream& os, vtkIndent indent) override;
40 
42 
45  int GetCellType() override { return VTK_CUBIC_LINE; }
46  int GetCellDimension() override { return 1; }
47  int GetNumberOfEdges() override { return 0; }
48  int GetNumberOfFaces() override { return 0; }
49  vtkCell* GetEdge(int) override { return nullptr; }
50  vtkCell* GetFace(int) override { return nullptr; }
51  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
52  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
53  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
54  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
55  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
56  double& dist2, double weights[]) override;
57  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
58  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
60  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
61  double* GetParametricCoords() override;
63 
68  double GetParametricDistance(const double pcoords[3]) override;
69 
74  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
75  vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
76  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
77 
81  int GetParametricCenter(double pcoords[3]) override;
82 
87  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
88  double pcoords[3], int& subId) override;
89 
90  static void InterpolationFunctions(const double pcoords[3], double weights[4]);
91  static void InterpolationDerivs(const double pcoords[3], double derivs[4]);
93 
97  void InterpolateFunctions(const double pcoords[3], double weights[4]) override
98  {
99  vtkCubicLine::InterpolationFunctions(pcoords, weights);
100  }
101  void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
102  {
103  vtkCubicLine::InterpolationDerivs(pcoords, derivs);
104  }
106 
107 protected:
109  ~vtkCubicLine() override;
110 
112  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
113 
114 private:
115  vtkCubicLine(const vtkCubicLine&) = delete;
116  void operator=(const vtkCubicLine&) = delete;
117 };
118 
119 //----------------------------------------------------------------------------
120 inline int vtkCubicLine::GetParametricCenter(double pcoords[3])
121 {
122 
123  pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
124  return 0;
125 }
126 
127 VTK_ABI_NAMESPACE_END
128 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:176
represent and manipulate cell attribute data
Definition: vtkCellData.h:31
abstract class to specify cell behavior
Definition: vtkCell.h:50
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
cell represents a cubic , isoparametric 1D line
Definition: vtkCubicLine.h:35
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
See the vtkCell API for descriptions of these methods.
static vtkCubicLine * New()
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:48
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
See the vtkCell API for descriptions of these methods.
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkCubicLine.h:97
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-line intersection.
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:47
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
int GetParametricCenter(double pcoords[3]) override
Return the center of the triangle in parametric coordinates.
Definition: vtkCubicLine.h:120
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *lines, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this line using scalar value provided.
~vtkCubicLine() override
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:45
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
See the vtkCell API for descriptions of these methods.
static void InterpolationDerivs(const double pcoords[3], double derivs[4])
vtkCell * GetFace(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:50
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkLine * Line
Definition: vtkCubicLine.h:111
void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkCubicLine.h:101
vtkCell * GetEdge(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:49
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:46
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
vtkDoubleArray * Scalars
Definition: vtkCubicLine.h:112
static void InterpolationFunctions(const double pcoords[3], double weights[4])
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
See the vtkCell API for descriptions of these methods.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:23
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:29
cell represents a 1D line
Definition: vtkLine.h:23
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:30
represent and manipulate 3D points
Definition: vtkPoints.h:29
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_CUBIC_LINE
Definition: vtkCellType.h:74
int vtkIdType
Definition: vtkType.h:315