ITK  4.4.0
Insight Segmentation and Registration Toolkit
Classes | Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
itk::fem::Element Class Referenceabstract

#include <itkFEMElementBase.h>

+ Inheritance diagram for itk::fem::Element:
+ Collaboration diagram for itk::fem::Element:

Detailed Description

Abstract base element class.

Derive this class to create new finite element classes. The storage of element parameters (geometry...) can't be implemented here, since we don't know yet, how much memory each element needs. Instead each derived class should take care of the memory management (declare appropriate data members) for the element parameters and provide access to these parameters (like nodes, materials...).

Derived classes must define the following class methods: GetIntegrationPointAndWeight GetNumberOfIntegrationPoints ShapeFunctions ShapeFunctionDerivatives GetLocalFromGlobalCoordinates JacobianDeterminant JacobianInverse PopulateEdgeIds

These are required for the loads to be properly applied properly to the element.

See Also
Element2DC0LinearLine
Element2DC0LinearQuadrilateral
Element2DC0LinearTriangular
Element2DC1Beam
Element2DC0QuadraticTriangular
Element3DC0LinearHexahedron
Element3DC0LinearTetrahedron
Element3DC0LinearTriangular
Element3DC0LinearTriangularLaplaceBeltrami

Definition at line 72 of file itkFEMElementBase.h.

Classes

class  Node
 

Public Types

enum  { InvalidDegreeOfFreedomID = 0xffffffff }
 
enum  { gaussMaxOrder = 10 }
 
typedef FEMPArray< ElementArrayType
 
typedef VectorContainer
< ElementIdentifier,
Element::Pointer
ArrayType1
 
typedef SmartPointer< const SelfConstPointer
 
typedef unsigned int DegreeOfFreedomIDType
 
typedef unsigned long ElementIdentifier
 
typedef double Float
 
typedef vnl_matrix< FloatMatrixType
 
typedef Node::ConstPointer NodeIDType
 
typedef SmartPointer< SelfPointer
 
typedef Element Self
 
typedef FEMLightObject Superclass
 
typedef vnl_vector< FloatVectorType
 
typedef FEMLightObject LoadType
 
typedef LoadType::Pointer LoadPointer
 
- Public Types inherited from itk::fem::FEMLightObject
typedef Self Baseclass
 
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef FEMLightObject Self
 
typedef itk::LightObject Superclass
 
- Public Types inherited from itk::LightObject
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef LightObject Self
 

Public Member Functions

virtual std::vector
< std::vector< int > > 
GetEdgeIds (void) const
 
virtual Float GetElementDeformationEnergy (MatrixType &LocalSolution) const
 
virtual VectorType GetGlobalFromLocalCoordinates (const VectorType &pt) const
 
virtual void GetIntegrationPointAndWeight (unsigned int i, VectorType &pt, Float &w, unsigned int order=0) const =0
 
virtual void GetLandmarkContributionMatrix (float eta, MatrixType &Le) const
 
virtual bool GetLocalFromGlobalCoordinates (const VectorType &globalPt, VectorType &localPt) const =0
 
virtual void GetMassMatrix (MatrixType &Me) const
 
virtual Material::ConstPointer GetMaterial (void) const
 
virtual void GetMaterialMatrix (MatrixType &D) const =0
 
virtual const char * GetNameOfClass () const
 
virtual NodeIDType GetNode (unsigned int n) const =0
 
virtual const VectorTypeGetNodeCoordinates (unsigned int n) const =0
 
virtual unsigned int GetNumberOfDegreesOfFreedom (void) const
 
virtual unsigned int GetNumberOfDegreesOfFreedomPerNode (void) const =0
 
virtual unsigned int GetNumberOfIntegrationPoints (unsigned int order=0) const =0
 
virtual unsigned int GetNumberOfNodes (void) const =0
 
virtual unsigned int GetNumberOfSpatialDimensions () const =0
 
virtual void GetStiffnessMatrix (MatrixType &Ke) const
 
virtual void GetStrainDisplacementMatrix (MatrixType &B, const MatrixType &shapeDgl) const =0
 
virtual VectorType GetStrainsAtPoint (const VectorType &pt, const Solution &sol, unsigned int index) const
 
virtual VectorType GetStressesAtPoint (const VectorType &pt, const VectorType &e, const Solution &sol, unsigned int index) const
 
virtual VectorType InterpolateSolution (const VectorType &pt, const Solution &sol, unsigned int solutionIndex=0) const
 
virtual Float InterpolateSolutionN (const VectorType &pt, const Solution &sol, unsigned int f, unsigned int solutionIndex=0) const
 
virtual void Jacobian (const VectorType &pt, MatrixType &J, const MatrixType *pshapeD=0) const
 
virtual Float JacobianDeterminant (const VectorType &pt, const MatrixType *pJ=0) const
 
virtual void JacobianInverse (const VectorType &pt, MatrixType &invJ, const MatrixType *pJ=0) const
 
virtual void PopulateEdgeIds (void)=0
 
virtual void SetMaterial (Material::ConstPointer)
 
virtual void SetNode (unsigned int n, NodeIDType node)=0
 
virtual void SetNode (unsigned int n, Node::Pointer node)
 
virtual void ShapeFunctionDerivatives (const VectorType &pt, MatrixType &shapeD) const =0
 
virtual void ShapeFunctionGlobalDerivatives (const VectorType &pt, MatrixType &shapeDgl, const MatrixType *pJ=0, const MatrixType *pshapeD=0) const
 
virtual VectorType ShapeFunctions (const VectorType &pt) const =0
 
DegreeOfFreedomIDType GetDegreeOfFreedom (unsigned int local_dof) const
 
- Public Member Functions inherited from itk::fem::FEMLightObject
int GetGlobalNumber () const
 
void SetGlobalNumber (int)
 
- Public Member Functions inherited from itk::LightObject
virtual Pointer CreateAnother () const
 
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
 itkCloneMacro (Self)
 
void Print (std::ostream &os, Indent indent=0) const
 
virtual void Register () const
 
virtual void SetReferenceCount (int)
 
virtual void UnRegister () const
 

Static Public Attributes

static const Float gaussPoint [gaussMaxOrder+1][gaussMaxOrder]
 
static const Float gaussWeight [gaussMaxOrder+1][gaussMaxOrder]
 

Protected Member Functions

virtual void PrintSelf (std::ostream &os, Indent indent) const
 
- Protected Member Functions inherited from itk::fem::FEMLightObject
 FEMLightObject ()
 
 FEMLightObject (const FEMLightObject &o)
 
virtual ~FEMLightObject ()
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Protected Attributes

std::vector< std::vector< int > > m_EdgeIds
 
- Protected Attributes inherited from itk::fem::FEMLightObject
int m_GlobalNumber
 
- Protected Attributes inherited from itk::LightObject
InternalReferenceCountType m_ReferenceCount
 
SimpleFastMutexLock m_ReferenceCountLock
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 
- Protected Types inherited from itk::LightObject
typedef int InternalReferenceCountType
 

Member Typedef Documentation

Array class that holds special pointers to the Element objects

Definition at line 94 of file itkFEMElementBase.h.

Definition at line 95 of file itkFEMElementBase.h.

Definition at line 79 of file itkFEMElementBase.h.

Type that stores global ID's of degrees of freedom.

Definition at line 126 of file itkFEMElementBase.h.

typedef unsigned long itk::fem::Element::ElementIdentifier

Definition at line 88 of file itkFEMElementBase.h.

typedef double itk::fem::Element::Float

Floating point type used in all Element classes.

Definition at line 82 of file itkFEMElementBase.h.

Easy and consistent access to LoadElement and LoadElement::Pointer type. This is a pointer to FEMLightObject to avoid cyclic references between LoadElement and Element classes. As a consequence whenever you need to use a pointer to LoadElement class within the element's declaration or definition, ALWAYS use this typedef instead. When calling the GetLoadVector(...) function from outside, you should ALWAYS first convert the argument to Element::LoadPointer. See code of function Solver::AssembleF(...) for more info.

Definition at line 119 of file itkFEMElementBase.h.

Easy and consistent access to LoadElement and LoadElement::Pointer type. This is a pointer to FEMLightObject to avoid cyclic references between LoadElement and Element classes. As a consequence whenever you need to use a pointer to LoadElement class within the element's declaration or definition, ALWAYS use this typedef instead. When calling the GetLoadVector(...) function from outside, you should ALWAYS first convert the argument to Element::LoadPointer. See code of function Solver::AssembleF(...) for more info.

Definition at line 118 of file itkFEMElementBase.h.

typedef vnl_matrix<Float> itk::fem::Element::MatrixType

Class used to store the element stiffness matrix

Definition at line 100 of file itkFEMElementBase.h.

Type that is used to store IDs of a node. It is a pointer to Node objects.

Definition at line 548 of file itkFEMElementBase.h.

Definition at line 78 of file itkFEMElementBase.h.

Standard class typedefs.

Definition at line 76 of file itkFEMElementBase.h.

Definition at line 77 of file itkFEMElementBase.h.

typedef vnl_vector<Float> itk::fem::Element::VectorType

Class to store the element load vector

Definition at line 105 of file itkFEMElementBase.h.

Member Enumeration Documentation

anonymous enum

Constant that represents an invalid DegreeOfFreedomID object. If a degree of freedom is assigned this value, this means that that no specific value was (yet) assigned to this DOF.

Enumerator
InvalidDegreeOfFreedomID 

Definition at line 133 of file itkFEMElementBase.h.

anonymous enum

Maximum supported order of 1D Gauss-Legendre integration. Integration points are defined for orders from 1 to gaussMaxOrder. Number of integration points is equal to the order of integration rule.

See Also
gaussPoint
Enumerator
gaussMaxOrder 

Definition at line 517 of file itkFEMElementBase.h.

Member Function Documentation

DegreeOfFreedomIDType itk::fem::Element::GetDegreeOfFreedom ( unsigned int  local_dof) const
inline

Convenient way to access IDs of degrees of freedom that are stored in node objects.

Parameters
local_dofLocal number of degree of freedom within an element.

Definition at line 420 of file itkFEMElementBase.h.

References GetNode(), GetNumberOfDegreesOfFreedom(), GetNumberOfDegreesOfFreedomPerNode(), and InvalidDegreeOfFreedomID.

virtual std::vector<std::vector<int> > itk::fem::Element::GetEdgeIds ( void  ) const
inlinevirtual

Access the edge ids vector. The vector in turn contains a list of edge ids.

Definition at line 706 of file itkFEMElementBase.h.

References m_EdgeIds.

virtual Float itk::fem::Element::GetElementDeformationEnergy ( MatrixType LocalSolution) const
virtual

Compute the physical energy, U, of the deformation (e.g. stress / strain ).

 T

U = u Ke u

The matrix LocalSolution contains the solution to use in the energy computation. Usually, this is the solution at the nodes.

virtual VectorType itk::fem::Element::GetGlobalFromLocalCoordinates ( const VectorType pt) const
virtual

Transforms the given local element coordinates into global.

Parameters
ptPoint in local element coordinates.
virtual void itk::fem::Element::GetIntegrationPointAndWeight ( unsigned int  i,
VectorType pt,
Float w,
unsigned int  order = 0 
) const
pure virtual

Methods related to numeric integration Computes the vector representing the i-th integration point in local element coordinates for a Gauss-Legendre numerical integration over the element domain. It also computes the weight at this integration point.

Optionally you can also specify the order of integration. If order is not specified, it defaults to 0, which means that the derived element should use the optimal integration order specific for that element.

Note
This function must be implemented in derived element classes, and is expected to provide valid integration points for up to gaussMaxOrder-th order of integration.
Parameters
iIntegration point number 0<=i<GetNumberOfIntegrationPoints()
ptReference to object of class VectorType that will hold the integration point.
wReference to Float variable that will hold the weight.
orderOrder of integration.
See Also
GetNumberOfIntegrationPoints()

Implemented in itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element3DC0LinearTriangular, and itk::fem::Element2DC0LinearLine.

virtual void itk::fem::Element::GetLandmarkContributionMatrix ( float  eta,
MatrixType Le 
) const
virtual

Compute and return landmark contribution to element stiffness matrix (Le) in global coordinate system.

b             T

int (1/eta)^2 N(x) N(x) dx a

where (eta ) is the landmark weight. Implementation is similar to GetMassMatrix.

virtual bool itk::fem::Element::GetLocalFromGlobalCoordinates ( const VectorType globalPt,
VectorType localPt 
) const
pure virtual

Transforms the given global element coordinates into local. Returns false if the point is outside.

Parameters
globalPtReference to vector containing a point in global (world) coordinates.
localPtReference to the vector that will store the local coordinate.

Implemented in itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element3DC0LinearTriangular, and itk::fem::Element2DC0LinearLine.

virtual void itk::fem::Element::GetMassMatrix ( MatrixType Me) const
virtual
virtual Material::ConstPointer itk::fem::Element::GetMaterial ( void  ) const
inlinevirtual

Return the pointer to the Material object used by the element. All derived classes, which use objects of Material class should override this method to provide access to the material from the base class.

Note
Derived Element classes don't have to use a material class, but since the majority of the final Element classes uses Material classes to specify phhysical constants that the element depends on, we provide this virtual function that enables easy access to this pointer from the base class. If the derived class does not override this function, the returned pointer is 0 by default, signaling that there is no Material object.
See Also
SetMaterial

Reimplemented in itk::fem::Element2DC1Beam, itk::fem::Element2DMembrane< Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< Element2DC0LinearTriangular >, itk::fem::Element1DStress< Element2DC0LinearLine >, itk::fem::Element2DStrain< Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTetrahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTriangular >, itk::fem::Element3DMembrane1DOF< Element3DC0LinearTriangular >, itk::fem::Element2DStress< Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< Element2DC0LinearTriangular >, itk::fem::Element3DStrain< Element3DC0LinearHexahedron >, and itk::fem::Element3DStrain< Element3DC0LinearTetrahedron >.

Definition at line 448 of file itkFEMElementBase.h.

virtual void itk::fem::Element::GetMaterialMatrix ( MatrixType D) const
pure virtual
virtual const char* itk::fem::Element::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::fem::FEMLightObject.

Reimplemented in itk::fem::Element2DC0LinearTriangularMembrane, itk::fem::Element2DC0QuadraticTriangularStrain, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangularStrain, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element2DC0QuadraticTriangularStress, itk::fem::Element2DC0LinearQuadrilateralMembrane, itk::fem::Element2DC0LinearQuadrilateralStress, itk::fem::Element2DC0LinearQuadrilateralStrain, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::ElementStd< 3, 2 >, itk::fem::ElementStd< 2, 2 >, itk::fem::ElementStd< 8, 3 >, itk::fem::ElementStd< 6, 2 >, itk::fem::ElementStd< 4, 3 >, itk::fem::ElementStd< 4, 2 >, itk::fem::ElementStd< 3, 3 >, itk::fem::Element2DMembrane< Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< Element2DC0LinearTriangular >, itk::fem::Element3DC0LinearTriangular, itk::fem::Element2DStrain< Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTetrahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTriangular >, itk::fem::Element3DC0LinearHexahedronMembrane, itk::fem::Element3DC0LinearHexahedronStrain, itk::fem::Element3DStrain< Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< Element3DC0LinearTetrahedron >, itk::fem::Element2DC0LinearTriangularStress, itk::fem::Element3DC0LinearTetrahedronMembrane, itk::fem::Element1DStress< Element2DC0LinearLine >, itk::fem::Element2DStress< Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< Element2DC0LinearTriangular >, itk::fem::Element3DMembrane1DOF< Element3DC0LinearTriangular >, itk::fem::Element3DC0LinearTetrahedronStrain, itk::fem::Element3DC0LinearTriangularLaplaceBeltrami, itk::fem::Element3DC0LinearTriangularMembrane, itk::fem::Element2DC1Beam, itk::fem::Element2DC0LinearLineStress, and itk::fem::Element2DC0LinearLine.

virtual NodeIDType itk::fem::Element::GetNode ( unsigned int  n) const
pure virtual
virtual const VectorType& itk::fem::Element::GetNodeCoordinates ( unsigned int  n) const
pure virtual

Return a vector of global coordinates of n-th node in an element.

Parameters
nLocal number of node. Must be 0 <= n < this->GetNumberOfNodes().

Implemented in itk::fem::ElementStd< 3, 2 >, itk::fem::ElementStd< 2, 2 >, itk::fem::ElementStd< 8, 3 >, itk::fem::ElementStd< 6, 2 >, itk::fem::ElementStd< 4, 3 >, itk::fem::ElementStd< 4, 2 >, and itk::fem::ElementStd< 3, 3 >.

virtual unsigned int itk::fem::Element::GetNumberOfDegreesOfFreedom ( void  ) const
inlinevirtual

Return the total number of degrees of freedom defined in a derived element class. By default this is equal to number of points in a cell multiplied by number of degrees of freedom at each point.

Definition at line 698 of file itkFEMElementBase.h.

References GetNumberOfDegreesOfFreedomPerNode(), and GetNumberOfNodes().

Referenced by GetDegreeOfFreedom().

virtual unsigned int itk::fem::Element::GetNumberOfDegreesOfFreedomPerNode ( void  ) const
pure virtual
virtual unsigned int itk::fem::Element::GetNumberOfIntegrationPoints ( unsigned int  order = 0) const
pure virtual

Returns total number of integration points, for given order of Gauss-Legendre numerical integration rule.

Note
This function must be implemented in derived element classes, and is expected to provide valid number of integration points for up to gaussMaxOrder-th order of integration.
See Also
GetIntegrationPointAndWeight()

Implemented in itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element3DC0LinearTriangular, and itk::fem::Element2DC0LinearLine.

virtual unsigned int itk::fem::Element::GetNumberOfNodes ( void  ) const
pure virtual
virtual unsigned int itk::fem::Element::GetNumberOfSpatialDimensions ( ) const
pure virtual

Returns the number of dimensions of space in which the element is defined. e.g. 2 for 2D elements, 3 for 3D... This is also equal to the size vector containing nodal coordinates.

Implemented in itk::fem::ElementStd< 3, 2 >, itk::fem::ElementStd< 2, 2 >, itk::fem::ElementStd< 8, 3 >, itk::fem::ElementStd< 6, 2 >, itk::fem::ElementStd< 4, 3 >, itk::fem::ElementStd< 4, 2 >, and itk::fem::ElementStd< 3, 3 >.

virtual void itk::fem::Element::GetStiffnessMatrix ( MatrixType Ke) const
virtual

Compute and return element stiffnes matrix (Ke) in global coordinate system. The base class provides a general implementation which only computes

b   T

int B(x) D B(x) dx a

using the Gaussian numeric integration method. The function calls GetIntegrationPointAndWeight() / GetNumberOfIntegrationPoints() to obtain the integration points. It also calls the GetStrainDisplacementMatrix() and GetMaterialMatrix() member functions.

Parameters
KeReference to the resulting stiffnes matrix.
Note
This is a very generic implementation of the stiffness matrix that is suitable for any problem/element definition. A specifc element may override this implementation with its own simple one.

Reimplemented in itk::fem::Element3DMembrane1DOF< Element3DC0LinearTriangular >, itk::fem::Element1DStress< Element2DC0LinearLine >, itk::fem::Element3DC0LinearTriangularLaplaceBeltrami, and itk::fem::Element2DC1Beam.

virtual void itk::fem::Element::GetStrainDisplacementMatrix ( MatrixType B,
const MatrixType shapeDgl 
) const
pure virtual
virtual VectorType itk::fem::Element::GetStrainsAtPoint ( const VectorType pt,
const Solution sol,
unsigned int  index 
) const
virtual
virtual VectorType itk::fem::Element::GetStressesAtPoint ( const VectorType pt,
const VectorType e,
const Solution sol,
unsigned int  index 
) const
virtual
virtual VectorType itk::fem::Element::InterpolateSolution ( const VectorType pt,
const Solution sol,
unsigned int  solutionIndex = 0 
) const
virtual

Return interpolated value of all unknown functions at given local point.

Parameters
ptPoint in local element coordinates.
solReference to the master solution object. This object is created by the Solver object when the whole FEM problem is solved and contains the values of unknown functions at nodes (degrees of freedom).
solutionIndexWe allow more than one solution vector to be stored - this selects which to use in interpolation.
virtual Float itk::fem::Element::InterpolateSolutionN ( const VectorType pt,
const Solution sol,
unsigned int  f,
unsigned int  solutionIndex = 0 
) const
virtual

Return interpolated value of f-th unknown function at given local point.

Parameters
ptPoint in local element coordinates.
solReference to the master solution object. This object is created by the Solver object when the whole FEM problem is solved and contains the values of unknown functions at nodes (degrees of freedom).
fNumber of unknown function to interpolate. Must be 0 <= f < GetNumberOfDegreesOfFreedomPerNode().
solutionIndexWe allow more than one solution vector to be stored - this selects which to use in interpolation.
virtual void itk::fem::Element::Jacobian ( const VectorType pt,
MatrixType J,
const MatrixType pshapeD = 0 
) const
virtual

Compute the Jacobian matrix of the transformation from local to global coordinates at a given local point.

A column in this matrix corresponds to a global coordinate, while a row corresponds to different local coordinates. E.g. element at row 2, col 3 contains derivative of the third global coordinate with respect to local coordinate number 2.

In order to compute the Jacobian, we normally need the shape function derivatives. If they are known, you should pass a pointer to an object of MatrixType that contains the shape function derivatives. If they are not known, pass null pointer and they will be computed automatically.

Parameters
ptPoint in local coordinates
Jreference to matrix object, which will contain the jacobian
pshapeDA pointer to derivatives of shape functions at point pt. If this pointer is 0, derivatives will be computed as necessary.

Reimplemented in itk::fem::Element2DC0LinearLine.

virtual Float itk::fem::Element::JacobianDeterminant ( const VectorType pt,
const MatrixType pJ = 0 
) const
virtual

Compute the determinant of the Jacobian matrix at a given point with respect to the local coordinate system.

Parameters
ptPoint in local element coordinates.
pJOptional pointer to Jacobian matrix computed at point pt. If this is set to 0, the Jacobian will be computed as necessary.

Reimplemented in itk::fem::Element2DC1Beam, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, and itk::fem::Element3DC0LinearTriangular.

virtual void itk::fem::Element::JacobianInverse ( const VectorType pt,
MatrixType invJ,
const MatrixType pJ = 0 
) const
virtual

Compute the inverse of the Jacobian matrix at a given point with respect to the local coordinate system.

Parameters
ptPoint in local element coordinates.
invJReference to the object of MatrixType that will store the computed inverse if Jacobian.
pJOptional pointer to Jacobian matrix computed at point pt. If this is set to 0, the Jacobian will be computed as necessary.

Reimplemented in itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, and itk::fem::Element3DC0LinearTriangular.

virtual void itk::fem::Element::PopulateEdgeIds ( void  )
pure virtual
virtual void itk::fem::Element::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
protectedvirtual

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 itk::fem::FEMLightObject.

Reimplemented in itk::fem::Element2DC1Beam, itk::fem::ElementStd< 3, 2 >, itk::fem::ElementStd< 2, 2 >, itk::fem::ElementStd< 8, 3 >, itk::fem::ElementStd< 6, 2 >, itk::fem::ElementStd< 4, 3 >, itk::fem::ElementStd< 4, 2 >, itk::fem::ElementStd< 3, 3 >, itk::fem::Element3DC0LinearTriangular, itk::fem::Element2DMembrane< Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< Element2DC0LinearTriangular >, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DStrain< Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< Element2DC0LinearTriangular >, itk::fem::Element1DStress< Element2DC0LinearLine >, itk::fem::Element3DMembrane< Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTetrahedron >, itk::fem::Element3DMembrane< Element3DC0LinearTriangular >, itk::fem::Element2DStress< Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< Element2DC0LinearTriangular >, itk::fem::Element3DMembrane1DOF< Element3DC0LinearTriangular >, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element3DStrain< Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< Element3DC0LinearTetrahedron >, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element2DC0LinearLine, itk::fem::Element2DC0QuadraticTriangularStrain, itk::fem::Element2DC0QuadraticTriangularStress, itk::fem::Element2DC0LinearQuadrilateralMembrane, itk::fem::Element2DC0LinearQuadrilateralStress, itk::fem::Element2DC0LinearTriangularMembrane, itk::fem::Element2DC0LinearQuadrilateralStrain, itk::fem::Element2DC0LinearTriangularStrain, itk::fem::Element3DC0LinearTriangularLaplaceBeltrami, itk::fem::Element2DC0LinearTriangularStress, itk::fem::Element3DC0LinearHexahedronMembrane, itk::fem::Element3DC0LinearHexahedronStrain, itk::fem::Element3DC0LinearTetrahedronMembrane, itk::fem::Element3DC0LinearTetrahedronStrain, itk::fem::Element2DC0LinearLineStress, and itk::fem::Element3DC0LinearTriangularMembrane.

virtual void itk::fem::Element::SetMaterial ( Material::ConstPointer  )
inlinevirtual
virtual void itk::fem::Element::SetNode ( unsigned int  n,
NodeIDType  node 
)
pure virtual
virtual void itk::fem::Element::SetNode ( unsigned int  n,
Node::Pointer  node 
)
inlinevirtual
virtual void itk::fem::Element::ShapeFunctionDerivatives ( const VectorType pt,
MatrixType shapeD 
) const
pure virtual

Compute the matrix of values of the shape functions derivatives with respect to local coordinates of this element at a given point.

A column in this matrix corresponds to a specific shape function, while a row corresponds to different local coordinates. E.g. element at row 2, col 3 contains derivative of shape function number 3 with respect to local coordinate number 2.

Parameters
ptPoint in local element coordinates.
shapeDReference to a matrix object, which will be filled with values of shape function derivatives.
See Also
ShapeFunctionGlobalDerivatives

Implemented in itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element3DC0LinearTriangular, and itk::fem::Element2DC0LinearLine.

virtual void itk::fem::Element::ShapeFunctionGlobalDerivatives ( const VectorType pt,
MatrixType shapeDgl,
const MatrixType pJ = 0,
const MatrixType pshapeD = 0 
) const
virtual

Compute matrix of shape function derivatives with respect to global coordinates.

A column in this matrix corresponds to a specific shape function, while a row corresponds to different global coordinates.

Parameters
ptPoint in local element coordinates.
shapeDglReference to a matrix object, which will be filled with values of shape function derivatives w.r.t. global (world) element coordinates.
pJOptional pointer to Jacobian matrix computed at point pt. If this is set to 0, the Jacobian will be computed as necessary.
pshapeDA pointer to derivatives of shape functions at point pt. If this pointer is 0, derivatives will be computed as necessary.
See Also
ShapeFunctionDerivatives
virtual VectorType itk::fem::Element::ShapeFunctions ( const VectorType pt) const
pure virtual

Returns a vector containing the values of all shape functions that define the geometry of a finite element at a given local point within an element.

Parameters
ptPoint in local element coordinates.

Implemented in itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element3DC0LinearTriangular, and itk::fem::Element2DC0LinearLine.

Member Data Documentation

const Float itk::fem::Element::gaussPoint[gaussMaxOrder+1][gaussMaxOrder]
static

Points for 1D Gauss-Legendre integration from -1 to 1. First index is order of integration, second index is the number of integration point.

Example: gaussPoint[4][2] returns third point of the 4th order integration rule. Subarray gaussPoint[0][...] does not provide useful information. It is there only to keep order index correct.

See Also
gaussWeight

Definition at line 530 of file itkFEMElementBase.h.

const Float itk::fem::Element::gaussWeight[gaussMaxOrder+1][gaussMaxOrder]
static

Weights for Gauss-Legendre integration.

See Also
gaussPoint

Definition at line 537 of file itkFEMElementBase.h.

std::vector<std::vector<int> > itk::fem::Element::m_EdgeIds
protected

Definition at line 726 of file itkFEMElementBase.h.

Referenced by GetEdgeIds().


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