Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itk::fem::Element Class Reference

#include <itkFEMElementBase.h>

Inheritance diagram for itk::fem::Element:

Inheritance graph
[legend]
Collaboration diagram for itk::fem::Element:

Collaboration graph
[legend]

List of all members.


Detailed Description

Abstract base element class.

Derive this class to create new finite element classes. All derived classes must define:

and optionally (if required):

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...).

Definition at line 75 of file itkFEMElementBase.h.


Public Types

enum  { InvalidDegreeOfFreedomID = 0xffffffff }
enum  { gaussMaxOrder = 10 }
typedef FEMPArray< ElementArrayType
typedef Self Baseclass
typedef const SelfConstPointer
typedef unsigned int DegreeOfFreedomIDType
typedef double Float
typedef vnl_matrix< FloatMatrixType
typedef Node::ConstPointer NodeIDType
typedef SelfPointer
typedef Element Self
typedef FEMLightObject Superclass
typedef vnl_vector< FloatVectorType
typedef LoadType::Pointer LoadPointer
typedef FEMLightObject LoadType

Public Member Functions

virtual int ClassID () const =0
virtual Baseclass::Pointer Clone () const =0
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 void GetLoadVector (LoadPointer l, VectorType &Fe) const =0
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 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 Read (std::istream &f, void *info)
virtual void SetMaterial (Material::ConstPointer)
virtual void SetNode (unsigned int n, NodeIDType node)=0
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
virtual void Write (std::ostream &f) const
DegreeOfFreedomIDType GetDegreeOfFreedom (unsigned int local_dof) const

Static Public Member Functions

static FEMLightObject::Pointer CreateFromStream (std::istream &f, void *info)
static void SkipWhiteSpace (std::istream &f)

Public Attributes

int GN

Static Public Attributes

static const Float gaussPoint [gaussMaxOrder+1][gaussMaxOrder]
static const Float gaussWeight [gaussMaxOrder+1][gaussMaxOrder]
static const std::string whitespaces

Classes

class  Node
 Class that stores information required to define a node. More...

Member Typedef Documentation

Array class that holds special pointers to the Element objects

Definition at line 88 of file itkFEMElementBase.h.

Store the base class typedef for easy access from derived classes. FEM_CLASS macro also expects this for the FEMOF...

Definition at line 64 of file itkFEMLightObject.h.

Const pointer or SmartPointer to an object.

Reimplemented from itk::fem::FEMLightObject.

Reimplemented in itk::fem::Element2DC0LinearLine, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element1DStress< itk::fem::Element2DC0LinearLine >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearTetrahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearTetrahedron >, 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 >, and itk::fem::ElementStd< 4, 2 >.

Definition at line 77 of file itkFEMElementBase.h.

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

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.

Reimplemented 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 >, and itk::fem::ElementStd< 4, 2 >.

Definition at line 111 of file itkFEMElementBase.h.

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

Reimplemented 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 >, and itk::fem::ElementStd< 4, 2 >.

Definition at line 505 of file itkFEMElementBase.h.

Pointer or SmartPointer to an object.

Reimplemented from itk::fem::FEMLightObject.

Reimplemented in itk::fem::Element2DC0LinearLine, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearHexahedron, itk::fem::Element3DC0LinearTetrahedron, itk::fem::Element1DStress< itk::fem::Element2DC0LinearLine >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearTetrahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearTetrahedron >, 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 >, and itk::fem::ElementStd< 4, 2 >.

Definition at line 77 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 125 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 473 of file itkFEMElementBase.h.


Member Function Documentation

virtual int itk::fem::FEMLightObject::ClassID (  )  const [pure virtual, inherited]

Returns the class ID of the object. This function is used to determine the class of the object without having to use the dynamic_cast operator.

Note:
Class must be registered with the FEMObjectFactory in order to create the class ID. Abstract classes don't define this function.

Implemented in itk::fem::FiniteDifferenceFunctionLoad< TMoving, TFixed >, itk::fem::Element2DC1Beam, itk::fem::Element::Node, itk::fem::ImageMetricLoad< TMoving, TFixed >, itk::fem::LoadBC, itk::fem::LoadBCMFC, itk::fem::LoadEdge, itk::fem::LoadElement, itk::fem::LoadGravConst, itk::fem::LoadLandmark, itk::fem::LoadNode, itk::fem::LoadPoint, itk::fem::LoadTest< TClass >, and itk::fem::MaterialLinearElasticity.

virtual Baseclass::Pointer itk::fem::FEMLightObject::Clone (  )  const [pure virtual, inherited]

static FEMLightObject::Pointer itk::fem::FEMLightObject::CreateFromStream ( std::istream &  f,
void *  info 
) [static, inherited]

Read object of any derived type from stream.

This static function creates an object of a class, which is derived from FEMLightObject. The class of object is first determined from the stream, then the object of that class is constructed using the FEMObjectFactory. Finally the data for this object is read from the stream, by calling the Read() member function.

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_dof Local number of degree of freedom within an element.

Definition at line 389 of file itkFEMElementBase.h.

References itk::fem::Element::Node::GetDegreeOfFreedom(), GetNode(), GetNumberOfDegreesOfFreedom(), GetNumberOfDegreesOfFreedomPerNode(), and InvalidDegreeOfFreedomID.

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:
pt Point 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:
i Integration point number 0<=i<GetNumberOfIntegrationPoints()
pt Reference to object of class VectorType that will hold the integration point.
w Reference to Float variable that will hold the weight.
order Order of integration.
See also:
GetNumberOfIntegrationPoints()

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 void itk::fem::Element::GetLoadVector ( LoadPointer  l,
VectorType Fe 
) const [pure virtual]

Compute and return the element load vector for a given external load. The class of load object determines the type of load acting on the elemnent. Basically this is the contribution of this element on the right side of the master matrix equation, due to the specified load. Returned vector includes only nodal forces that correspond to the given Load object.

Visitor design pattern is used in the loads implementation. This function only selects and calls the proper function based on the given class of load object. The code that performs the actual conversion to the corresponding nodal loads is defined elswhere.

Note:
Each derived class must implement its own version of this function. This is automated by calling the LOAD_FUNCTION() macro within the class declaration (in the public: block).
For example on how to define specific element load, see funtion LoadImplementationPoint_Bar2D.

Note:
: Before a load can be applied to an element, the function that implements a load must be registered with the VisitorDispactcher class.
Parameters:
l Pointer to a load object.
Fe Reference to vector object that will store nodal forces.
See also:
VisitorDispatcher

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:
globalPt Reference to vector containing a point in global (world) coordinates.
localPt Reference to the vector that will store the local coordinate.

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

virtual void itk::fem::Element::GetMassMatrix ( MatrixType Me  )  const [virtual]

virtual Material::ConstPointer itk::fem::Element::GetMaterial ( void   )  const [inline, virtual]

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::Element1DStress< itk::fem::Element2DC0LinearLine >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearTetrahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearHexahedron >, and itk::fem::Element3DStrain< itk::fem::Element3DC0LinearTetrahedron >.

Definition at line 412 of file itkFEMElementBase.h.

virtual void itk::fem::Element::GetMaterialMatrix ( MatrixType D  )  const [pure virtual]

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:
n Local 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 >, and itk::fem::ElementStd< 4, 2 >.

virtual unsigned int itk::fem::Element::GetNumberOfDegreesOfFreedom ( void   )  const [inline, virtual]

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 650 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::Element2DC0LinearLine, itk::fem::Element2DC0LinearQuadrilateral, itk::fem::Element2DC0LinearTriangular, itk::fem::Element2DC0QuadraticTriangular, itk::fem::Element2DC1Beam, itk::fem::Element3DC0LinearHexahedron, and itk::fem::Element3DC0LinearTetrahedron.

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 >, and itk::fem::ElementStd< 4, 2 >.

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:
Ke Reference 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::Element2DC1Beam, and itk::fem::Element1DStress< itk::fem::Element2DC0LinearLine >.

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:
pt Point in local element coordinates.
sol Reference 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).
solutionIndex We 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:
pt Point in local element coordinates.
sol Reference 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).
f Number of unknown function to interpolate. Must be 0 <= f < GetNumberOfDegreesOfFreedomPerNode().
solutionIndex We 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:
pt Point in local coordinates
J referece to matrix object, which will contain the jacobian
pshapeD A 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:
pt Point in local element coordinates.
pJ Optional 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::Element2DC1Beam.

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:
pt Point in local element coordinates.
invJ Reference to the object of MatrixType that will store the computed inverse if Jacobian.
pJ Optional 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, and itk::fem::Element2DC0QuadraticTriangular.

virtual void itk::fem::FEMLightObject::Read ( std::istream &  f,
void *  info 
) [virtual, inherited]

Read an object data from input stream. Call this member to initialize the data members in the current object by reading data from provided input stream. Derived classes should first call the the parent's read function, to initialize the data from parent. Note that you must manually create the object of desired type using the FEMObjectFactory before you can call read function (this is pretty obvious). In this class only the global number is read from file. Derived classes may require some additional info in order to perform the reading. Pack this info in an object and pass a pointer to it in the info parameter. If you need runtime typechecking, use a polymorphic class and dynamic_cast operator inside the implementation of Read.

Reimplemented in itk::fem::Element2DC1Beam, itk::fem::Element::Node, itk::fem::LoadBC, itk::fem::LoadBCMFC, itk::fem::LoadEdge, itk::fem::LoadElement, itk::fem::LoadGravConst, itk::fem::LoadLandmark, itk::fem::LoadNode, itk::fem::LoadTest< TClass >, itk::fem::MaterialLinearElasticity, itk::fem::Element1DStress< itk::fem::Element2DC0LinearLine >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearTetrahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearTetrahedron >, 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 >, and itk::fem::ElementStd< 4, 2 >.

virtual void itk::fem::Element::SetMaterial ( Material::ConstPointer   )  [inline, virtual]

virtual void itk::fem::Element::SetNode ( unsigned int  n,
NodeIDType  node 
) [pure virtual]

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:
pt Point in local element coordinates.
shapeD Reference to a matrix object, which will be filled with values of shape function derivatives.
See also:
ShapeFunctionGlobalDerivatives

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

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:
pt Point in local element coordinates.
shapeDgl Reference to a matrix object, which will be filled with values of shape function derivatives w.r.t. global (world) element coordinates.
pJ Optional pointer to Jacobian matrix computed at point pt. If this is set to 0, the Jacobian will be computed as necessary.
pshapeD A 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:
pt Point in local element coordinates.

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

static void itk::fem::FEMLightObject::SkipWhiteSpace ( std::istream &  f  )  [static, inherited]

Helper function that skips all the whitespace and comments in an input stream.

virtual void itk::fem::FEMLightObject::Write ( std::ostream &  f  )  const [virtual, inherited]

Write an object to the output stream. Call this member to write the data members in the current object to the output stream. Here we also need to know which derived class we actually are, so that we can write the class name. The class name is obtained by calling the virtual ClassID() member function and passing the result to the FEMObjectFactory.

Implementations of Write member funtion in derived classes should first call the parent's implementation of Write and finaly write whatever they need.

Reimplemented in itk::fem::Element2DC1Beam, itk::fem::Element::Node, itk::fem::LoadBC, itk::fem::LoadBCMFC, itk::fem::LoadEdge, itk::fem::LoadElement, itk::fem::LoadGravConst, itk::fem::LoadLandmark, itk::fem::LoadNode, itk::fem::LoadTest< TClass >, itk::fem::MaterialLinearElasticity, itk::fem::Element1DStress< itk::fem::Element2DC0LinearLine >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DMembrane< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStrain< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStrain< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearQuadrilateral >, itk::fem::Element2DStress< itk::fem::Element2DC0QuadraticTriangular >, itk::fem::Element2DStress< itk::fem::Element2DC0LinearTriangular >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DMembrane< itk::fem::Element3DC0LinearTetrahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearHexahedron >, itk::fem::Element3DStrain< itk::fem::Element3DC0LinearTetrahedron >, 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 >, and itk::fem::ElementStd< 4, 2 >.


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 486 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 493 of file itkFEMElementBase.h.

Global number of an object (ID of an object) In general the ID's are required to be unique only within a specific type of derived classes (Elements, Nodes, ...) If the GN is not required, it can be ignored. (normally you need the GN when writing or reading objects to/from stream.

Definition at line 165 of file itkFEMLightObject.h.

Referenced by itk::fem::FEMLightObject::FEMLightObject().

const std::string itk::fem::FEMLightObject::whitespaces [static, inherited]

Const string of all whitespace characters. This string is used by SkipWhiteSpace function.

Definition at line 135 of file itkFEMLightObject.h.


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

Generated at Thu May 28 18:44:04 2009 for ITK by doxygen 1.5.5 written by Dimitri van Heesch, © 1997-2000