VTK
vtkLagrangeCurve.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLagrangeCurve.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 =========================================================================*/
15 // .NAME vtkLagrangeCurve
16 // .SECTION Description
17 // .SECTION See Also
18 
19 #ifndef vtkLagrangeCurve_h
20 #define vtkLagrangeCurve_h
21 
22 #include "vtkCommonDataModelModule.h" // For export macro
23 #include "vtkNonLinearCell.h"
24 #include "vtkSmartPointer.h" // For member variable.
25 #include "vtkCellType.h" // For GetCellType.
26 #include "vtkNew.h" // For member variable.
27 
28 class vtkCellData;
29 class vtkDoubleArray;
30 class vtkIdList;
31 class vtkLine;
32 class vtkPointData;
33 class vtkPoints;
34 class vtkVector3d;
35 class vtkVector3i;
36 
37 class VTKCOMMONDATAMODEL_EXPORT vtkLagrangeCurve : public vtkNonLinearCell
38 {
39 public:
40  static vtkLagrangeCurve* New();
42  void PrintSelf(ostream& os, vtkIndent indent) override;
43 
44  int GetCellType() override { return VTK_LAGRANGE_CURVE; }
45  int GetCellDimension() override { return 1; }
46  int RequiresInitialization() override { return 1; }
47  int GetNumberOfEdges() override { return 0; }
48  int GetNumberOfFaces() override { return 0; }
49  vtkCell* GetEdge(int) override { return nullptr; }
50  vtkCell* GetFace(int) override { return nullptr; }
51 
52 
53  void Initialize() override;
54 
55  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
56  int EvaluatePosition(const double x[3], double closestPoint[3],
57  int& subId, double pcoords[3],
58  double& dist2, double weights[]) override;
59  void EvaluateLocation(
60  int& subId, const double pcoords[3], double x[3],
61  double* weights) override;
62  void Contour(
63  double value, vtkDataArray* cellScalars,
65  vtkCellArray* lines, vtkCellArray* polys,
66  vtkPointData* inPd, vtkPointData* outPd,
67  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
68  void Clip(
69  double value, vtkDataArray* cellScalars,
71  vtkPointData* inPd, vtkPointData* outPd,
72  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd,
73  int insideOut) override;
74  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
75  double x[3], double pcoords[3], int& subId) override;
76  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
77  void Derivatives(
78  int subId, const double pcoords[3], const double* values,
79  int dim, double* derivs) override;
80  double* GetParametricCoords() override;
81  int GetParametricCenter(double center[3]) override;
82 
83  double GetParametricDistance(const double pcoords[3]) override;
84 
85  const int* GetOrder();
86  int GetOrder(int i) { return this->GetOrder()[i]; }
87 
88  void InterpolateFunctions(const double pcoords[3], double* weights) override;
89  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
90 
91  bool SubCellCoordinatesFromId(vtkVector3i& ijk, int subId);
92  bool SubCellCoordinatesFromId(int& i, int subId);
93  int PointIndexFromIJK(int i, int, int);
94  bool TransformApproxToCellParams(int subCell, double* pcoords);
95 
96 protected:
98  ~vtkLagrangeCurve() override;
99 
100  vtkLine* GetApprox();
101  void PrepareApproxData(vtkPointData* pd, vtkCellData* cd, vtkIdType cellId, vtkDataArray* cellScalars);
102  vtkLine* GetApproximateLine(
103  int subId, vtkDataArray* scalarsIn = nullptr, vtkDataArray* scalarsOut = nullptr);
104 
105  int Order[2];
114 
115 private:
116  vtkLagrangeCurve(const vtkLagrangeCurve&) = delete;
117  void operator=(const vtkLagrangeCurve&) = delete;
118 };
119 
121 {
122  center[0] = 0.5;
123  center[1] = center[2] = 0.0;
124  return 0;
125 }
126 
127 #endif // vtkLagrangeCurve_h
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
vtkSmartPointer< vtkPoints > PointParametricCoordinates
represent and manipulate point attribute data
Definition: vtkPointData.h:37
int GetNumberOfFaces() override
Return the number of faces in the cell.
int GetNumberOfEdges() override
Return the number of edges in the cell.
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
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.
virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this leve...
Definition: vtkCell.h:357
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:347
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkNew< vtkDoubleArray > Scalars
vtkNew< vtkPoints > TmpPts
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...
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:360
vtkNew< vtkIdList > TmpIds
dynamic, self-adjusting array of double
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:35
abstract class to specify cell behavior
Definition: vtkCell.h:59
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:36
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
int GetCellType() override
Return the type of cell.
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.
vtkSmartPointer< vtkPointData > ApproxPD
vtkSmartPointer< vtkCellData > ApproxCD
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:50
int GetParametricCenter(double center[3]) override
Return center of the cell in parametric coordinates.
vtkNew< vtkDoubleArray > CellScalars
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.
virtual void Initialize()
Definition: vtkCell.h:114
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.
int RequiresInitialization() override
Some cells require initialization prior to access.
represent and manipulate 3D points
Definition: vtkPoints.h:39
vtkSmartPointer< vtkLine > Approx