VTK
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkLagrangeTriangle Class Reference

A 2D cell that represents an arbitrary order Lagrange triangle. More...

#include <vtkLagrangeTriangle.h>

Inheritance diagram for vtkLagrangeTriangle:
[legend]
Collaboration diagram for vtkLagrangeTriangle:
[legend]

Public Types

typedef vtkNonLinearCell Superclass
 
- Public Types inherited from vtkNonLinearCell
typedef vtkCell Superclass
 
- Public Types inherited from vtkCell
typedef vtkObject Superclass
 

Public Member Functions

virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkLagrangeTriangleNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int GetCellType () override
 Return the type of cell. More...
 
int GetCellDimension () override
 Return the topological dimensional of the cell (0,1,2, or 3). More...
 
int RequiresInitialization () override
 Some cells require initialization prior to access. More...
 
int GetNumberOfEdges () override
 Return the number of edges in the cell. More...
 
int GetNumberOfFaces () override
 Return the number of faces in the cell. More...
 
vtkCellGetEdge (int edgeId) override
 Return the edge cell from the edgeId of the cell. More...
 
vtkCellGetFace (int) override
 Return the face cell from the faceId of the cell. More...
 
void Initialize () override
 
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 inside or outside of the cell. More...
 
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; evaluate parametric coordinates, sub-cell id (!=0 only if cell is composite), distance squared of point x[3] to cell (in particular, the sub-cell indicated), closest point on cell to x[3] (unless closestPoint is null, in which case, the closest point and dist2 are not found), and interpolation weights in cell. More...
 
void EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 Determine global coordinate (x[3]) from subId and parametric coordinates. More...
 
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. More...
 
void Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
 Cut (or clip) the cell based on the input cellScalars and the specified value. More...
 
int IntersectWithLine (const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
 Intersect with a ray. More...
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 Generate simplices of proper dimension. More...
 
void JacobianInverse (const double pcoords[3], double **inverse, double *derivs)
 
void Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
 Compute derivatives given cell subId and parametric coordinates. More...
 
doubleGetParametricCoords () override
 Return a contiguous array of parametric coordinates of the points defining this cell. More...
 
int GetParametricCenter (double pcoords[3]) override
 Return center of the cell in parametric coordinates. More...
 
double GetParametricDistance (const double pcoords[3]) override
 Return the distance of the parametric coordinate provided to the cell. More...
 
void InterpolateFunctions (const double pcoords[3], double *weights) override
 
void InterpolateDerivs (const double pcoords[3], double *derivs) override
 
vtkIdType GetOrder () const
 
vtkIdType ComputeOrder ()
 
void ToBarycentricIndex (vtkIdType index, vtkIdType *bindex)
 
vtkIdType ToIndex (const vtkIdType *bindex)
 
- Public Member Functions inherited from vtkNonLinearCell
vtkNonLinearCellNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int IsLinear () override
 Non-linear cells require special treatment (tessellation) when converting to graphics primitives (during mapping). More...
 
- Public Member Functions inherited from vtkCell
vtkCellNewInstance () const
 
void Initialize (int npts, vtkIdType *pts, vtkPoints *p)
 Initialize cell from outside with point ids and point coordinates specified. More...
 
void Initialize (int npts, vtkPoints *p)
 Initialize the cell with point coordinates specified. More...
 
virtual void ShallowCopy (vtkCell *c)
 Copy this cell by reference counting the internal data structures. More...
 
virtual void DeepCopy (vtkCell *c)
 Copy this cell by completely copying internal data structures. More...
 
virtual int IsExplicitCell ()
 Explicit cells require additional representational information beyond the usual cell type and connectivity list information. More...
 
virtual int RequiresExplicitFaceRepresentation ()
 Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods). More...
 
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
 
virtual vtkIdTypeGetFaces ()
 
vtkPointsGetPoints ()
 Get the point coordinates for the cell. More...
 
vtkIdType GetNumberOfPoints ()
 Return the number of points in the cell. More...
 
vtkIdListGetPointIds ()
 Return the list of point ids defining the cell. More...
 
vtkIdType GetPointId (int ptId)
 For cell point i, return the actual point id. More...
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
doubleGetBounds ()
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
double GetLength2 ()
 Compute Length squared of cell (i.e., bounding box diagonal squared). More...
 
virtual int IsPrimaryCell ()
 Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet). More...
 
virtual void InterpolateFunctions (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. More...
 
virtual void InterpolateDerivs (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on. More...
 
virtual void DebugOff ()
 Turn debugging output off. More...
 
bool GetDebug ()
 Get the value of the debug flag. More...
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag. More...
 
virtual void Modified ()
 Update the modification time for this object. More...
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time. More...
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
vtkTypeBool HasObserver (unsigned long event)
 
vtkTypeBool HasObserver (const char *event)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Allow user to set the AbortFlagOn() with the return value of the callback method. More...
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
virtual void Delete ()
 Delete a VTK object. More...
 
virtual void FastDelete ()
 Delete a reference to this object. More...
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream. More...
 
virtual void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object). More...
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object). More...
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
void PrintRevisions (ostream &)
 Legacy. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 

Static Public Member Functions

static vtkLagrangeTriangleNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkLagrangeTriangleSafeDownCast (vtkObjectBase *o)
 
static int MaximumOrder ()
 
static int MaximumNumberOfPoints ()
 
static void ComputeParametricCoords (double *, vtkIdType)
 
static void BarycentricIndex (vtkIdType index, vtkIdType *bindex, vtkIdType order)
 
static vtkIdType Index (const vtkIdType *bindex, vtkIdType order)
 
static double eta (vtkIdType n, vtkIdType chi, double sigma)
 
static double d_eta (vtkIdType n, vtkIdType chi, double sigma)
 
- Static Public Member Functions inherited from vtkNonLinearCell
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkNonLinearCellSafeDownCast (vtkObjectBase *o)
 
- Static Public Member Functions inherited from vtkCell
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkCellSafeDownCast (vtkObjectBase *o)
 
- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes. More...
 
static void SetGlobalWarningDisplay (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static int GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 Return 1 if this class type is the same type of (or a subclass of) the named class. More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkLagrangeTriangle ()
 
 ~vtkLagrangeTriangle () override
 
vtkIdType GetNumberOfSubtriangles () const
 
vtkIdType ComputeNumberOfSubtriangles ()
 
void SubtriangleBarycentricPointIndices (vtkIdType cellIndex, vtkIdType(&pointBIndices)[3][3])
 
- Protected Member Functions inherited from vtkNonLinearCell
 vtkNonLinearCell ()
 
 ~vtkNonLinearCell () override
 
- Protected Member Functions inherited from vtkCell
 vtkCell ()
 
 ~vtkCell () override
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
 ~vtkObject () override
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr)
 These methods allow a command to exclusively grab all events. More...
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events. More...
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkLagrangeCurveEdge
 
vtkTriangleFace
 
vtkDoubleArrayScalars
 
vtkIdType Order
 
vtkIdType NumberOfSubtriangles
 
doubleParametricCoordinates
 
vtkIdType EdgeIds [VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1]
 
vtkIdType BarycentricIndexMap [3 *MAX_POINTS]
 
vtkIdType IndexMap [(VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1)*(VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1)]
 
vtkIdType SubtriangleIndexMap [9 *MAX_SUBTRIANGLES]
 
- Protected Attributes inherited from vtkCell
double Bounds [6]
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkCell
vtkPointsPoints
 
vtkIdListPointIds
 

Detailed Description

A 2D cell that represents an arbitrary order Lagrange triangle.

vtkLagrangeTriangle is a concrete implementation of vtkCell to represent a 2D triangle using Lagrange shape functions of user specified order.

The number of points in a Lagrange cell determines the order over which they are iterated relative to the parametric coordinate system of the cell. The first points that are reported are vertices. They appear in the same order in which they would appear in linear cells. Mid-edge points are reported next. They are reported in sequence. For two- and three-dimensional (3D) cells, the following set of points to be reported are face points. Finally, 3D cells report points interior to their volume.

Tests:
vtkLagrangeTriangle (Tests)

Definition at line 52 of file vtkLagrangeTriangle.h.

Member Typedef Documentation

Definition at line 56 of file vtkLagrangeTriangle.h.

Constructor & Destructor Documentation

vtkLagrangeTriangle::vtkLagrangeTriangle ( )
protected
vtkLagrangeTriangle::~vtkLagrangeTriangle ( )
overrideprotected

Member Function Documentation

static vtkLagrangeTriangle* vtkLagrangeTriangle::New ( )
static
static vtkTypeBool vtkLagrangeTriangle::IsTypeOf ( const char *  type)
static
virtual vtkTypeBool vtkLagrangeTriangle::IsA ( const char *  name)
virtual

Return 1 if this class is the same type of (or a subclass of) the named class.

Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkNonLinearCell.

static vtkLagrangeTriangle* vtkLagrangeTriangle::SafeDownCast ( vtkObjectBase o)
static
virtual vtkObjectBase* vtkLagrangeTriangle::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkNonLinearCell.

vtkLagrangeTriangle* vtkLagrangeTriangle::NewInstance ( ) const
void vtkLagrangeTriangle::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

Methods invoked by print to print information about the object including superclasses.

Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from vtkCell.

int vtkLagrangeTriangle::GetCellType ( )
inlineoverridevirtual

Return the type of cell.

Implements vtkCell.

Definition at line 59 of file vtkLagrangeTriangle.h.

int vtkLagrangeTriangle::GetCellDimension ( )
inlineoverridevirtual

Return the topological dimensional of the cell (0,1,2, or 3).

Implements vtkCell.

Definition at line 60 of file vtkLagrangeTriangle.h.

int vtkLagrangeTriangle::RequiresInitialization ( )
inlineoverridevirtual

Some cells require initialization prior to access.

For example, they may have to triangulate themselves or set up internal data structures.

Reimplemented from vtkCell.

Definition at line 61 of file vtkLagrangeTriangle.h.

int vtkLagrangeTriangle::GetNumberOfEdges ( )
inlineoverridevirtual

Return the number of edges in the cell.

Implements vtkCell.

Definition at line 62 of file vtkLagrangeTriangle.h.

int vtkLagrangeTriangle::GetNumberOfFaces ( )
inlineoverridevirtual

Return the number of faces in the cell.

Implements vtkCell.

Definition at line 63 of file vtkLagrangeTriangle.h.

vtkCell* vtkLagrangeTriangle::GetEdge ( int  edgeId)
overridevirtual

Return the edge cell from the edgeId of the cell.

Implements vtkCell.

vtkCell* vtkLagrangeTriangle::GetFace ( int  faceId)
inlineoverridevirtual

Return the face cell from the faceId of the cell.

Implements vtkCell.

Definition at line 65 of file vtkLagrangeTriangle.h.

void vtkLagrangeTriangle::Initialize ( )
overridevirtual

Reimplemented from vtkCell.

static int vtkLagrangeTriangle::MaximumOrder ( )
inlinestatic

Definition at line 69 of file vtkLagrangeTriangle.h.

static int vtkLagrangeTriangle::MaximumNumberOfPoints ( )
inlinestatic

Definition at line 70 of file vtkLagrangeTriangle.h.

int vtkLagrangeTriangle::CellBoundary ( int  subId,
const double  pcoords[3],
vtkIdList pts 
)
overridevirtual

Given parametric coordinates of a point, return the closest cell boundary, and whether the point is inside or outside of the cell.

The cell boundary is defined by a list of points (pts) that specify a face (3D cell), edge (2D cell), or vertex (1D cell). If the return value of the method is != 0, then the point is inside the cell.

Implements vtkCell.

int vtkLagrangeTriangle::EvaluatePosition ( const double  x[3],
double  closestPoint[3],
int subId,
double  pcoords[3],
double dist2,
double  weights[] 
)
overridevirtual

Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; evaluate parametric coordinates, sub-cell id (!=0 only if cell is composite), distance squared of point x[3] to cell (in particular, the sub-cell indicated), closest point on cell to x[3] (unless closestPoint is null, in which case, the closest point and dist2 are not found), and interpolation weights in cell.

(The number of weights is equal to the number of points defining the cell). Note: on rare occasions a -1 is returned from the method. This means that numerical error has occurred and all data returned from this method should be ignored. Also, inside/outside is determine parametrically. That is, a point is inside if it satisfies parametric limits. This can cause problems for cells of topological dimension 2 or less, since a point in 3D can project onto the cell within parametric limits but be "far" from the cell. Thus the value dist2 may be checked to determine true in/out.

Implements vtkCell.

void vtkLagrangeTriangle::EvaluateLocation ( int subId,
const double  pcoords[3],
double  x[3],
double weights 
)
overridevirtual

Determine global coordinate (x[3]) from subId and parametric coordinates.

Also returns interpolation weights. (The number of weights is equal to the number of points in the cell.)

Implements vtkCell.

void vtkLagrangeTriangle::Contour ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray verts,
vtkCellArray lines,
vtkCellArray polys,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd 
)
overridevirtual

Generate contouring primitives.

The scalar list cellScalars are scalar values at each cell point. The point locator is essentially a points list that merges points as they are inserted (i.e., prevents duplicates). Contouring primitives can be vertices, lines, or polygons. It is possible to interpolate point data along the edge by providing input and output point data - if outPd is nullptr, then no interpolation is performed. Also, if the output cell data is non-nullptr, the cell data from the contoured cell is passed to the generated contouring primitives. (Note: the CopyAllocate() method must be invoked on both the output cell and point data. The cellId refers to the cell from which the cell data is copied.)

Implements vtkCell.

void vtkLagrangeTriangle::Clip ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray connectivity,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd,
int  insideOut 
)
overridevirtual

Cut (or clip) the cell based on the input cellScalars and the specified value.

The output of the clip operation will be one or more cells of the same topological dimension as the original cell. The flag insideOut controls what part of the cell is considered inside - normally cell points whose scalar value is greater than "value" are considered inside. If insideOut is on, this is reversed. Also, if the output cell data is non-nullptr, the cell data from the clipped cell is passed to the generated contouring primitives. (Note: the CopyAllocate() method must be invoked on both the output cell and point data. The cellId refers to the cell from which the cell data is copied.)

Implements vtkCell.

int vtkLagrangeTriangle::IntersectWithLine ( const double  p1[3],
const double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId 
)
overridevirtual

Intersect with a ray.

Return parametric coordinates (both line and cell) and global intersection coordinates, given ray definition p1[3], p2[3] and tolerance tol. The method returns non-zero value if intersection occurs. A parametric distance t between 0 and 1 along the ray representing the intersection point, the point coordinates x[3] in data coordinates and also pcoords[3] in parametric coordinates. subId is the index within the cell if a composed cell like a triangle strip.

Implements vtkCell.

int vtkLagrangeTriangle::Triangulate ( int  index,
vtkIdList ptIds,
vtkPoints pts 
)
overridevirtual

Generate simplices of proper dimension.

If cell is 3D, tetrahedron are generated; if 2D triangles; if 1D lines; if 0D points. The form of the output is a sequence of points, each n+1 points (where n is topological cell dimension) defining a simplex. The index is a parameter that controls which triangulation to use (if more than one is possible). If numerical degeneracy encountered, 0 is returned, otherwise 1 is returned. This method does not insert new points: all the points that define the simplices are the points that define the cell.

Implements vtkCell.

void vtkLagrangeTriangle::JacobianInverse ( const double  pcoords[3],
double **  inverse,
double derivs 
)
void vtkLagrangeTriangle::Derivatives ( int  subId,
const double  pcoords[3],
const double values,
int  dim,
double derivs 
)
overridevirtual

Compute derivatives given cell subId and parametric coordinates.

The values array is a series of data value(s) at the cell points. There is a one-to-one correspondence between cell point and data value(s). Dim is the number of data values per cell point. Derivs are derivatives in the x-y-z coordinate directions for each data value. Thus, if computing derivatives for a scalar function in a hexahedron, dim=1, 8 values are supplied, and 3 deriv values are returned (i.e., derivatives in x-y-z directions). On the other hand, if computing derivatives of velocity (vx,vy,vz) dim=3, 24 values are supplied ((vx,vy,vz)1, (vx,vy,vz)2, ....()8), and 9 deriv values are returned ((d(vx)/dx),(d(vx)/dy),(d(vx)/dz), (d(vy)/dx),(d(vy)/dy), (d(vy)/dz), (d(vz)/dx),(d(vz)/dy),(d(vz)/dz)).

Implements vtkCell.

double* vtkLagrangeTriangle::GetParametricCoords ( )
overridevirtual

Return a contiguous array of parametric coordinates of the points defining this cell.

In other words, (px,py,pz, px,py,pz, etc..) The coordinates are ordered consistent with the definition of the point ordering for the cell. This method returns a non-nullptr pointer when the cell is a primary type (i.e., IsPrimaryCell() is true). Note that 3D parametric coordinates are returned no matter what the topological dimension of the cell.

Reimplemented from vtkCell.

static void vtkLagrangeTriangle::ComputeParametricCoords ( double ,
vtkIdType   
)
static
int vtkLagrangeTriangle::GetParametricCenter ( double  pcoords[3])
overridevirtual

Return center of the cell in parametric coordinates.

Note that the parametric center is not always located at (0.5,0.5,0.5). The return value is the subId that the center is in (if a composite cell). If you want the center in x-y-z space, invoke the EvaluateLocation() method.

Reimplemented from vtkCell.

double vtkLagrangeTriangle::GetParametricDistance ( const double  pcoords[3])
overridevirtual

Return the distance of the parametric coordinate provided to the cell.

If inside the cell, a distance of zero is returned. This is used during picking to get the correct cell picked. (The tolerance will occasionally allow cells to be picked who are not really intersected "inside" the cell.)

Reimplemented from vtkCell.

void vtkLagrangeTriangle::InterpolateFunctions ( const double  pcoords[3],
double weights 
)
override
void vtkLagrangeTriangle::InterpolateDerivs ( const double  pcoords[3],
double derivs 
)
override
vtkIdType vtkLagrangeTriangle::GetOrder ( ) const
inline

Definition at line 107 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::ComputeOrder ( )
void vtkLagrangeTriangle::ToBarycentricIndex ( vtkIdType  index,
vtkIdType bindex 
)
vtkIdType vtkLagrangeTriangle::ToIndex ( const vtkIdType bindex)
static void vtkLagrangeTriangle::BarycentricIndex ( vtkIdType  index,
vtkIdType bindex,
vtkIdType  order 
)
static
static vtkIdType vtkLagrangeTriangle::Index ( const vtkIdType bindex,
vtkIdType  order 
)
static
static double vtkLagrangeTriangle::eta ( vtkIdType  n,
vtkIdType  chi,
double  sigma 
)
static
static double vtkLagrangeTriangle::d_eta ( vtkIdType  n,
vtkIdType  chi,
double  sigma 
)
static
vtkIdType vtkLagrangeTriangle::GetNumberOfSubtriangles ( ) const
inlineprotected

Definition at line 124 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::ComputeNumberOfSubtriangles ( )
protected
void vtkLagrangeTriangle::SubtriangleBarycentricPointIndices ( vtkIdType  cellIndex,
vtkIdType(&)  pointBIndices[3][3] 
)
protected

Member Data Documentation

vtkLagrangeCurve* vtkLagrangeTriangle::Edge
protected

Definition at line 133 of file vtkLagrangeTriangle.h.

vtkTriangle* vtkLagrangeTriangle::Face
protected

Definition at line 134 of file vtkLagrangeTriangle.h.

vtkDoubleArray* vtkLagrangeTriangle::Scalars
protected

Definition at line 135 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::Order
protected

Definition at line 136 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::NumberOfSubtriangles
protected

Definition at line 137 of file vtkLagrangeTriangle.h.

double* vtkLagrangeTriangle::ParametricCoordinates
protected

Definition at line 138 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::EdgeIds[VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1]
protected

Definition at line 140 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::BarycentricIndexMap[3 *MAX_POINTS]
protected

Definition at line 141 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::IndexMap[(VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1)*(VTK_LAGRANGE_TRIANGLE_MAX_ORDER+1)]
protected

Definition at line 143 of file vtkLagrangeTriangle.h.

vtkIdType vtkLagrangeTriangle::SubtriangleIndexMap[9 *MAX_SUBTRIANGLES]
protected

Definition at line 144 of file vtkLagrangeTriangle.h.


The documentation for this class was generated from the following file: