VTK  9.2.6
vtkCubicLine.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCubicLine.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 =========================================================================*/
36 #ifndef vtkCubicLine_h
37 #define vtkCubicLine_h
38 
39 #include "vtkCommonDataModelModule.h" // For export macro
40 #include "vtkNonLinearCell.h"
41 
42 class vtkLine;
43 class vtkDoubleArray;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkCubicLine : public vtkNonLinearCell
46 {
47 public:
48  static vtkCubicLine* New();
50  void PrintSelf(ostream& os, vtkIndent indent) override;
51 
53 
56  int GetCellType() override { return VTK_CUBIC_LINE; }
57  int GetCellDimension() override { return 1; }
58  int GetNumberOfEdges() override { return 0; }
59  int GetNumberOfFaces() override { return 0; }
60  vtkCell* GetEdge(int) override { return nullptr; }
61  vtkCell* GetFace(int) override { return nullptr; }
62  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
63  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
64  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
65  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
66  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
67  double& dist2, double weights[]) override;
68  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
69  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
70  void Derivatives(
71  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
72  double* GetParametricCoords() override;
74 
79  double GetParametricDistance(const double pcoords[3]) override;
80 
85  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
86  vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
87  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
88 
92  int GetParametricCenter(double pcoords[3]) override;
93 
98  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
99  double pcoords[3], int& subId) override;
100 
101  static void InterpolationFunctions(const double pcoords[3], double weights[4]);
102  static void InterpolationDerivs(const double pcoords[3], double derivs[4]);
104 
108  void InterpolateFunctions(const double pcoords[3], double weights[4]) override
109  {
110  vtkCubicLine::InterpolationFunctions(pcoords, weights);
111  }
112  void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
113  {
114  vtkCubicLine::InterpolationDerivs(pcoords, derivs);
115  }
117 
118 protected:
119  vtkCubicLine();
120  ~vtkCubicLine() override;
121 
123  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
124 
125 private:
126  vtkCubicLine(const vtkCubicLine&) = delete;
127  void operator=(const vtkCubicLine&) = delete;
128 };
129 
130 //----------------------------------------------------------------------------
131 inline int vtkCubicLine::GetParametricCenter(double pcoords[3])
132 {
133 
134  pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
135  return 0;
136 }
137 
138 #endif
vtkLine * Line
Definition: vtkCubicLine.h:122
void InterpolateDerivs(const double pcoords[3], double derivs[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
Definition: vtkCubicLine.h:112
represent and manipulate point attribute data
Definition: vtkPointData.h:41
vtkCell * GetFace(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:61
represent and manipulate cell attribute data
Definition: vtkCellData.h:41
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
abstract superclass for non-linear cells
static void InterpolationDerivs(const double pcoords[3], double derivs[4])
int vtkIdType
Definition: vtkType.h:332
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:58
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:56
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
dynamic, self-adjusting array of double
cell represents a cubic , isoparametric 1D line
Definition: vtkCubicLine.h:45
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
cell represents a 1D line
Definition: vtkLine.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:60
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
a simple class to control print indentation
Definition: vtkIndent.h:39
static void InterpolationFunctions(const double pcoords[3], double weights[4])
vtkDoubleArray * Scalars
Definition: vtkCubicLine.h:123
list of point or cell ids
Definition: vtkIdList.h:33
int GetParametricCenter(double pcoords[3]) override
Return the center of the triangle in parametric coordinates.
Definition: vtkCubicLine.h:131
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
Definition: vtkCubicLine.h:108
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:59
object to represent cell connectivity
Definition: vtkCellArray.h:186
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
vtkCell * GetEdge(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:60
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkCubicLine.h:57
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
represent and manipulate 3D points
Definition: vtkPoints.h:39