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

cell represents a parabolic, isoparametric edge More...

#include <vtkQuadraticEdge.h>

Inheritance diagram for vtkQuadraticEdge:
[legend]
Collaboration diagram for vtkQuadraticEdge:
[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...
 
vtkQuadraticEdgeNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int GetCellType () override
 Implement the vtkCell API. More...
 
int GetCellDimension () override
 Return the topological dimensional of the cell (0,1,2, or 3). 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) 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...
 
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...
 
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...
 
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...
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 Generate simplices of proper dimension. More...
 
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...
 
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. More...
 
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. More...
 
int GetParametricCenter (double pcoords[3]) override
 Return the center of the quadratic tetra in parametric coordinates. More...
 
void InterpolateFunctions (const double pcoords[3], double weights[3]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
void InterpolateDerivs (const double pcoords[3], double derivs[3]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
- Public Member Functions inherited from vtkNonLinearCell
vtkNonLinearCellNewInstance () const
 
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 RequiresInitialization ()
 Some cells require initialization prior to access. More...
 
virtual void Initialize ()
 
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 double GetParametricDistance (const double pcoords[3])
 Return the distance of the parametric coordinate provided to the cell. 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 vtkQuadraticEdgeNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkQuadraticEdgeSafeDownCast (vtkObjectBase *o)
 
static void InterpolationFunctions (const double pcoords[3], double weights[3])
 
static void InterpolationDerivs (const double pcoords[3], double derivs[3])
 
- 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
 
 vtkQuadraticEdge ()
 
 ~vtkQuadraticEdge () override
 
- 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

vtkLineLine
 
vtkDoubleArrayScalars
 
- 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

cell represents a parabolic, isoparametric edge

vtkQuadraticEdge is a concrete implementation of vtkNonLinearCell to represent a one-dimensional, 3-nodes, isoparametric parabolic line. The interpolation is the standard finite element, quadratic isoparametric shape function. The cell includes a mid-edge node. The ordering of the three points defining the cell is point ids (0,1,2) where id #2 is the midedge node.

See also
vtkQuadraticTriangle vtkQuadraticTetra vtkQuadraticWedge vtkQuadraticQuad vtkQuadraticHexahedron vtkQuadraticPyramid
Tests:
vtkQuadraticEdge (Tests)

Definition at line 43 of file vtkQuadraticEdge.h.

Member Typedef Documentation

Definition at line 47 of file vtkQuadraticEdge.h.

Constructor & Destructor Documentation

vtkQuadraticEdge::vtkQuadraticEdge ( )
protected
vtkQuadraticEdge::~vtkQuadraticEdge ( )
overrideprotected

Member Function Documentation

static vtkQuadraticEdge* vtkQuadraticEdge::New ( )
static
static vtkTypeBool vtkQuadraticEdge::IsTypeOf ( const char *  type)
static
virtual vtkTypeBool vtkQuadraticEdge::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 vtkQuadraticEdge* vtkQuadraticEdge::SafeDownCast ( vtkObjectBase o)
static
virtual vtkObjectBase* vtkQuadraticEdge::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkNonLinearCell.

vtkQuadraticEdge* vtkQuadraticEdge::NewInstance ( ) const
void vtkQuadraticEdge::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 vtkNonLinearCell.

int vtkQuadraticEdge::GetCellType ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 54 of file vtkQuadraticEdge.h.

int vtkQuadraticEdge::GetCellDimension ( )
inlineoverridevirtual

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

Implements vtkCell.

Definition at line 55 of file vtkQuadraticEdge.h.

int vtkQuadraticEdge::GetNumberOfEdges ( )
inlineoverridevirtual

Return the number of edges in the cell.

Implements vtkCell.

Definition at line 56 of file vtkQuadraticEdge.h.

int vtkQuadraticEdge::GetNumberOfFaces ( )
inlineoverridevirtual

Return the number of faces in the cell.

Implements vtkCell.

Definition at line 57 of file vtkQuadraticEdge.h.

vtkCell* vtkQuadraticEdge::GetEdge ( int  edgeId)
inlineoverridevirtual

Return the edge cell from the edgeId of the cell.

Implements vtkCell.

Definition at line 58 of file vtkQuadraticEdge.h.

vtkCell* vtkQuadraticEdge::GetFace ( int  faceId)
inlineoverridevirtual

Return the face cell from the faceId of the cell.

Implements vtkCell.

Definition at line 59 of file vtkQuadraticEdge.h.

int vtkQuadraticEdge::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.

void vtkQuadraticEdge::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.

int vtkQuadraticEdge::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 vtkQuadraticEdge::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.

int vtkQuadraticEdge::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 vtkQuadraticEdge::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* vtkQuadraticEdge::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.

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

Clip this edge using scalar value provided.

Like contouring, except that it cuts the edge to produce linear line segments.

Implements vtkCell.

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

Line-edge intersection.

Intersection has to occur within [0,1] parametric coordinates and with specified tolerance.

Implements vtkCell.

int vtkQuadraticEdge::GetParametricCenter ( double  pcoords[3])
inlineoverridevirtual

Return the center of the quadratic tetra in parametric coordinates.

Reimplemented from vtkCell.

Definition at line 134 of file vtkQuadraticEdge.h.

static void vtkQuadraticEdge::InterpolationFunctions ( const double  pcoords[3],
double  weights[3] 
)
static
static void vtkQuadraticEdge::InterpolationDerivs ( const double  pcoords[3],
double  derivs[3] 
)
static
void vtkQuadraticEdge::InterpolateFunctions ( const double  pcoords[3],
double  weights[3] 
)
inlineoverride

Compute the interpolation functions/derivatives (aka shape functions/derivatives)

Definition at line 112 of file vtkQuadraticEdge.h.

void vtkQuadraticEdge::InterpolateDerivs ( const double  pcoords[3],
double  derivs[3] 
)
inlineoverride

Compute the interpolation functions/derivatives (aka shape functions/derivatives)

Definition at line 116 of file vtkQuadraticEdge.h.

Member Data Documentation

vtkLine* vtkQuadraticEdge::Line
protected

Definition at line 126 of file vtkQuadraticEdge.h.

vtkDoubleArray* vtkQuadraticEdge::Scalars
protected

Definition at line 127 of file vtkQuadraticEdge.h.


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