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
33 #ifndef vtkCubicLine_h
34 #define vtkCubicLine_h
35 
36 #include "vtkCommonDataModelModule.h" // For export macro
37 #include "vtkNonLinearCell.h"
38 
39 VTK_ABI_NAMESPACE_BEGIN
40 class vtkLine;
41 class vtkDoubleArray;
42 
43 class VTKCOMMONDATAMODEL_EXPORT vtkCubicLine : public vtkNonLinearCell
44 {
45 public:
46  static vtkCubicLine* New();
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
51 
54  int GetCellType() override { return VTK_CUBIC_LINE; }
55  int GetCellDimension() override { return 1; }
56  int GetNumberOfEdges() override { return 0; }
57  int GetNumberOfFaces() override { return 0; }
58  vtkCell* GetEdge(int) override { return nullptr; }
59  vtkCell* GetFace(int) override { return nullptr; }
60  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
61  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
62  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
63  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
64  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
65  double& dist2, double weights[]) override;
66  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
67  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
69  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
70  double* GetParametricCoords() override;
72 
77  double GetParametricDistance(const double pcoords[3]) override;
78 
83  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
84  vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
85  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
86 
90  int GetParametricCenter(double pcoords[3]) override;
91 
96  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
97  double pcoords[3], int& subId) override;
98 
99  static void InterpolationFunctions(const double pcoords[3], double weights[4]);
100  static void InterpolationDerivs(const double pcoords[3], double derivs[4]);
102 
106  void InterpolateFunctions(const double pcoords[3], double weights[4]) override
107  {
108  vtkCubicLine::InterpolationFunctions(pcoords, weights);
109  }
110  void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
111  {
112  vtkCubicLine::InterpolationDerivs(pcoords, derivs);
113  }
115 
116 protected:
118  ~vtkCubicLine() override;
119 
121  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
122 
123 private:
124  vtkCubicLine(const vtkCubicLine&) = delete;
125  void operator=(const vtkCubicLine&) = delete;
126 };
127 
128 //----------------------------------------------------------------------------
129 inline int vtkCubicLine::GetParametricCenter(double pcoords[3])
130 {
131 
132  pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
133  return 0;
134 }
135 
136 VTK_ABI_NAMESPACE_END
137 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:185
represent and manipulate cell attribute data
Definition: vtkCellData.h:40
abstract class to specify cell behavior
Definition: vtkCell.h:59
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
cell represents a cubic , isoparametric 1D line
Definition: vtkCubicLine.h:44
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:57
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:106
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:56
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:129
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:54
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:59
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkLine * Line
Definition: vtkCubicLine.h:120
void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkCubicLine.h:110
vtkCell * GetEdge(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:58
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:55
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
vtkDoubleArray * Scalars
Definition: vtkCubicLine.h:121
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:54
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:38
cell represents a 1D line
Definition: vtkLine.h:32
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:39
represent and manipulate 3D points
Definition: vtkPoints.h:38
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_CUBIC_LINE
Definition: vtkCellType.h:83
int vtkIdType
Definition: vtkType.h:315