VTK  9.2.6
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 =========================================================================*/
34 #ifndef vtkQuadraticEdge_h
35 #define vtkQuadraticEdge_h
36 
37 #include "vtkCommonDataModelModule.h" // For export macro
38 #include "vtkNonLinearCell.h"
39 
40 class vtkLine;
41 class vtkDoubleArray;
42 
43 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticEdge : public vtkNonLinearCell
44 {
45 public:
46  static vtkQuadraticEdge* New();
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
54  int GetCellType() override { return VTK_QUADRATIC_EDGE; }
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 
61  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
62  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
63  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
64  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
65  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
66  double& dist2, double weights[]) override;
67  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
68  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
69  void Derivatives(
70  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
71  double* GetParametricCoords() override;
72 
77  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
78  vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
79  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
80 
85  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
86  double pcoords[3], int& subId) override;
87 
91  int GetParametricCenter(double pcoords[3]) override;
92 
93  static void InterpolationFunctions(const double pcoords[3], double weights[3]);
94  static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
96 
100  void InterpolateFunctions(const double pcoords[3], double weights[3]) override
101  {
103  }
104  void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
105  {
106  vtkQuadraticEdge::InterpolationDerivs(pcoords, derivs);
107  }
109 
110 protected:
112  ~vtkQuadraticEdge() override;
113 
115  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
116 
117 private:
118  vtkQuadraticEdge(const vtkQuadraticEdge&) = delete;
119  void operator=(const vtkQuadraticEdge&) = delete;
120 };
121 //----------------------------------------------------------------------------
122 inline int vtkQuadraticEdge::GetParametricCenter(double pcoords[3])
123 {
124  pcoords[0] = 0.5;
125  pcoords[1] = pcoords[2] = 0.;
126  return 0;
127 }
128 
129 #endif
vtkDoubleArray * Scalars
represent and manipulate point attribute data
Definition: vtkPointData.h:41
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
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.
int GetNumberOfFaces() override
Return the number of faces in the cell.
int GetCellType() override
Implement the vtkCell API.
abstract superclass for non-linear cells
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.
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
int GetNumberOfEdges() override
Return the number of edges in 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
list of point or cell ids
Definition: vtkIdList.h:33
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
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.
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...
object to represent cell connectivity
Definition: vtkCellArray.h:186
static void InterpolationFunctions(const double pcoords[3], double weights[3])
cell represents a parabolic, isoparametric edge
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.
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
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.
static void InterpolationDerivs(const double pcoords[3], double derivs[3])
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.
void InterpolateFunctions(const double pcoords[3], double weights[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate 3D points
Definition: vtkPoints.h:39
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic tetra in parametric coordinates.