VTK  9.1.0
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=========================================================================*/
40#ifndef vtkQuadraticEdge_h
41#define vtkQuadraticEdge_h
42
43#include "vtkCommonDataModelModule.h" // For export macro
44#include "vtkNonLinearCell.h"
45
46class vtkLine;
47class vtkDoubleArray;
48
49class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticEdge : public vtkNonLinearCell
50{
51public:
54 void PrintSelf(ostream& os, vtkIndent indent) override;
55
60 int GetCellType() override { return VTK_QUADRATIC_EDGE; }
61 int GetCellDimension() override { return 1; }
62 int GetNumberOfEdges() override { return 0; }
63 int GetNumberOfFaces() override { return 0; }
64 vtkCell* GetEdge(int) override { return nullptr; }
65 vtkCell* GetFace(int) override { return nullptr; }
66
67 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
68 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
69 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
70 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
71 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
72 double& dist2, double weights[]) override;
73 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
74 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
76 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
77 double* GetParametricCoords() 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
91 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
92 double pcoords[3], int& subId) override;
93
97 int GetParametricCenter(double pcoords[3]) override;
98
99 static void InterpolationFunctions(const double pcoords[3], double weights[3]);
100 static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
102
106 void InterpolateFunctions(const double pcoords[3], double weights[3]) override
107 {
109 }
110 void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
111 {
113 }
115
116protected:
119
121 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
122
123private:
124 vtkQuadraticEdge(const vtkQuadraticEdge&) = delete;
125 void operator=(const vtkQuadraticEdge&) = delete;
126};
127//----------------------------------------------------------------------------
128inline int vtkQuadraticEdge::GetParametricCenter(double pcoords[3])
129{
130 pcoords[0] = 0.5;
131 pcoords[1] = pcoords[2] = 0.;
132 return 0;
133}
134
135#endif
object to represent cell connectivity
Definition: vtkCellArray.h:190
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:67
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:40
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:43
cell represents a 1D line
Definition: vtkLine.h:40
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:42
represent and manipulate 3D points
Definition: vtkPoints.h:43
cell represents a parabolic, isoparametric edge
int GetCellType() override
Implement the vtkCell API.
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)
static vtkQuadraticEdge * New()
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)
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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.
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
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
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
int GetNumberOfEdges() override
Return the number of edges in the cell.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
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:65
int vtkIdType
Definition: vtkType.h:332