VTK  9.2.5
vtkQuadraticEdge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticEdge.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 =========================================================================*/
53 #ifndef vtkQuadraticEdge_h
54 #define vtkQuadraticEdge_h
55 
56 #include "vtkCommonDataModelModule.h" // For export macro
57 #include "vtkNonLinearCell.h"
58 
59 class vtkLine;
60 class vtkDoubleArray;
61 
62 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticEdge : public vtkNonLinearCell
63 {
64 public:
65  static vtkQuadraticEdge* New();
67  void PrintSelf(ostream& os, vtkIndent indent) override;
68 
73  int GetCellType() override { return VTK_QUADRATIC_EDGE; }
74  int GetCellDimension() override { return 1; }
75  int GetNumberOfEdges() override { return 0; }
76  int GetNumberOfFaces() override { return 0; }
77  vtkCell* GetEdge(int) override { return nullptr; }
78  vtkCell* GetFace(int) override { return nullptr; }
79 
80  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
81  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
82  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
83  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
84  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
85  double& dist2, double weights[]) override;
86  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
87  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
89  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
90  double* GetParametricCoords() override;
91 
96  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
97  vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
98  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
99 
104  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
105  double pcoords[3], int& subId) override;
106 
110  int GetParametricCenter(double pcoords[3]) override;
111 
112  static void InterpolationFunctions(const double pcoords[3], double weights[3]);
113  static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
115 
119  void InterpolateFunctions(const double pcoords[3], double weights[3]) override
120  {
122  }
123  void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
124  {
125  vtkQuadraticEdge::InterpolationDerivs(pcoords, derivs);
126  }
128 
129 protected:
131  ~vtkQuadraticEdge() override;
132 
134  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
135 
136 private:
137  vtkQuadraticEdge(const vtkQuadraticEdge&) = delete;
138  void operator=(const vtkQuadraticEdge&) = delete;
139 };
140 //----------------------------------------------------------------------------
141 inline int vtkQuadraticEdge::GetParametricCenter(double pcoords[3])
142 {
143  pcoords[0] = 0.5;
144  pcoords[1] = pcoords[2] = 0.;
145  return 0;
146 }
147 
148 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:143
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:119
cell represents a 1D line
Definition: vtkLine.h:143
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
cell represents a parabolic, isoparametric edge
int GetCellType() override
Implement the vtkCell API.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
static void InterpolationFunctions(const double pcoords[3], double weights[3])
void InterpolateFunctions(const double pcoords[3], double weights[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkDoubleArray * Scalars
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
static void InterpolationDerivs(const double pcoords[3], double derivs[3])
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic tetra in parametric coordinates.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
~vtkQuadraticEdge() override
int GetNumberOfEdges() override
Return the number of edges in the cell.
static vtkQuadraticEdge * New()
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
int GetNumberOfFaces() override
Return the number of faces in the cell.
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 edge using scalar value provided.
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
Generate contouring primitives.
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_EDGE
Definition: vtkCellType.h:104
int vtkIdType
Definition: vtkType.h:332