MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Types | Public Member Functions | Protected Attributes | List of all members
mfem::FiniteElement Class Referenceabstract

Abstract class for Finite Elements. More...

#include <fe.hpp>

Inheritance diagram for mfem::FiniteElement:
[legend]
Collaboration diagram for mfem::FiniteElement:
[legend]

Public Types

enum  { SCALAR, VECTOR }
 Enumeration for RangeType. More...
 
enum  { VALUE, INTEGRAL, H_DIV, H_CURL }
 

Public Member Functions

 FiniteElement (int D, int G, int Do, int O, int F=FunctionSpace::Pk)
 
int GetDim () const
 Returns the space dimension for the finite element. More...
 
int GetGeomType () const
 Returns the geometry type: More...
 
int GetDof () const
 Returns the degrees of freedom in the FE space. More...
 
int GetOrder () const
 Returns the order of the finite element. More...
 
int Space () const
 Returns the type of space on each element. More...
 
int GetRangeType () const
 
int GetMapType () const
 
virtual void CalcShape (const IntegrationPoint &ip, Vector &shape) const =0
 
virtual void CalcDShape (const IntegrationPoint &ip, DenseMatrix &dshape) const =0
 
const IntegrationRuleGetNodes () const
 
virtual void CalcVShape (const IntegrationPoint &ip, DenseMatrix &shape) const
 
virtual void CalcVShape (ElementTransformation &Trans, DenseMatrix &shape) const
 
virtual void CalcDivShape (const IntegrationPoint &ip, Vector &divshape) const
 
virtual void CalcCurlShape (const IntegrationPoint &ip, DenseMatrix &curl_shape) const
 
virtual void GetFaceDofs (int face, int **dofs, int *ndofs) const
 
virtual void CalcHessian (const IntegrationPoint &ip, DenseMatrix &h) const
 
virtual void GetLocalInterpolation (ElementTransformation &Trans, DenseMatrix &I) const
 
virtual void Project (Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
 
virtual void Project (VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
 
virtual void ProjectDelta (int vertex, Vector &dofs) const
 
virtual void Project (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
 
virtual void ProjectGrad (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
 
virtual void ProjectCurl (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
 
virtual void ProjectDiv (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
 
virtual ~FiniteElement ()
 

Protected Attributes

int Dim
 
int GeomType
 
int Dof
 
int Order
 
int FuncSpace
 
int RangeType
 
int MapType
 
IntegrationRule Nodes
 

Detailed Description

Abstract class for Finite Elements.

Definition at line 44 of file fe.hpp.

Member Enumeration Documentation

anonymous enum

Enumeration for RangeType.

Enumerator
SCALAR 
VECTOR 

Definition at line 52 of file fe.hpp.

anonymous enum

Enumeration for MapType: defines how reference functions, uh(xh), are mapped to functions on a general element, u(x):

VALUE u(x) = uh(xh) INTEGRAL u(x) = (1/w) * uh(xh) H_DIV u(x) = (J/w) * uh(xh) H_CURL u(x) = J^{-t} * uh(xh) (square J) H_CURL u(x) = J*(J^t*J)^{-1} * uh(xh) (general J)

where

x = T(xh) is the image of the reference point xh ("x hat"), J = J(xh) is the Jacobian matrix of the transformation T, and w = w(xh) = / det(J), for square J, . \ det(J^t*J)^{1/2}, for general J, is the transformation weight factor.

Enumerator
VALUE 
INTEGRAL 
H_DIV 
H_CURL 

Definition at line 71 of file fe.hpp.

Constructor & Destructor Documentation

mfem::FiniteElement::FiniteElement ( int  D,
int  G,
int  Do,
int  O,
int  F = FunctionSpace::Pk 
)

Construct Finite Element with given (D)im - space dimension, (G)eomType - geometry type (of type Geometry::Type), (Do)f - degrees of freedom in the FE space, and (O)rder - order of the FE space (F)uncSpace- type of space on each element

Definition at line 22 of file fe.cpp.

virtual mfem::FiniteElement::~FiniteElement ( )
inlinevirtual

Definition at line 196 of file fe.hpp.

Member Function Documentation

void mfem::FiniteElement::CalcCurlShape ( const IntegrationPoint ip,
DenseMatrix curl_shape 
) const
virtual

pure virtual function which evaluates the values of the curl all shape functions at a given point ip and stores them in the matrix curl_shape (Dof x Dim) so that each row contains the curl of one shape function

Reimplemented in mfem::ND_TriangleElement, mfem::ND_TetrahedronElement, mfem::ND_QuadrilateralElement, mfem::ND_HexahedronElement, mfem::Nedelec1TetFiniteElement, and mfem::Nedelec1HexFiniteElement.

Definition at line 51 of file fe.cpp.

void mfem::FiniteElement::CalcDivShape ( const IntegrationPoint ip,
Vector divshape 
) const
virtual
virtual void mfem::FiniteElement::CalcDShape ( const IntegrationPoint ip,
DenseMatrix dshape 
) const
pure virtual

pure virtual function which evaluates the values of all partial derivatives of all shape functions at a given point ip and stores them in the matrix dshape (Dof x Dim) so that each row contains the derivatives of one shape function

Implemented in mfem::NURBS3DFiniteElement, mfem::NURBS2DFiniteElement, mfem::NURBS1DFiniteElement, mfem::L2Pos_TetrahedronElement, mfem::L2_TetrahedronElement, mfem::L2Pos_TriangleElement, mfem::L2_TriangleElement, mfem::L2Pos_HexahedronElement, mfem::L2_HexahedronElement, mfem::L2Pos_QuadrilateralElement, mfem::L2_QuadrilateralElement, mfem::L2Pos_SegmentElement, mfem::L2_SegmentElement, mfem::H1_TetrahedronElement, mfem::H1_TriangleElement, mfem::H1Pos_HexahedronElement, mfem::H1Pos_QuadrilateralElement, mfem::H1Pos_SegmentElement, mfem::H1_HexahedronElement, mfem::H1_QuadrilateralElement, mfem::H1_SegmentElement, mfem::RotTriLinearHexFiniteElement, mfem::RefinedTriLinear3DFiniteElement, mfem::RefinedBiLinear2DFiniteElement, mfem::RefinedLinear3DFiniteElement, mfem::RefinedLinear2DFiniteElement, mfem::RefinedLinear1DFiniteElement, mfem::LagrangeHexFiniteElement, mfem::P0HexFiniteElement, mfem::P0TetFiniteElement, mfem::P1TetNonConfFiniteElement, mfem::Lagrange1DFiniteElement, mfem::P2SegmentFiniteElement, mfem::P1SegmentFiniteElement, mfem::P0SegmentFiniteElement, mfem::CrouzeixRaviartQuadFiniteElement, mfem::CrouzeixRaviartFiniteElement, mfem::TriLinear3DFiniteElement, mfem::Quadratic3DFiniteElement, mfem::Linear3DFiniteElement, mfem::P0QuadFiniteElement, mfem::P0TriangleFiniteElement, mfem::Cubic3DFiniteElement, mfem::Cubic2DFiniteElement, mfem::Cubic1DFiniteElement, mfem::BiCubic2DFiniteElement, mfem::GaussBiQuad2DFiniteElement, mfem::BiQuadPos2DFiniteElement, mfem::BiQuad2DFiniteElement, mfem::GaussQuad2DFiniteElement, mfem::Quad2DFiniteElement, mfem::QuadPos1DFiniteElement, mfem::Quad1DFiniteElement, mfem::P1OnQuadFiniteElement, mfem::GaussBiLinear2DFiniteElement, mfem::GaussLinear2DFiniteElement, mfem::BiLinear2DFiniteElement, mfem::Linear2DFiniteElement, mfem::Linear1DFiniteElement, and mfem::PointFiniteElement.

void mfem::FiniteElement::CalcHessian ( const IntegrationPoint ip,
DenseMatrix h 
) const
virtual

each row of h contains the upper triangular part of the hessian of one shape function; the order in 2D is {u_xx, u_xy, u_yy}

Reimplemented in mfem::Cubic2DFiniteElement, mfem::BiCubic2DFiniteElement, mfem::Quad2DFiniteElement, and mfem::BiLinear2DFiniteElement.

Definition at line 63 of file fe.cpp.

virtual void mfem::FiniteElement::CalcShape ( const IntegrationPoint ip,
Vector shape 
) const
pure virtual

pure virtual function which evaluates the values of all shape functions at a given point ip and stores them in the vector shape of dimension Dof

Implemented in mfem::NURBS3DFiniteElement, mfem::NURBS2DFiniteElement, mfem::NURBS1DFiniteElement, mfem::ND_SegmentElement, mfem::L2Pos_TetrahedronElement, mfem::L2_TetrahedronElement, mfem::L2Pos_TriangleElement, mfem::L2_TriangleElement, mfem::L2Pos_HexahedronElement, mfem::L2_HexahedronElement, mfem::L2Pos_QuadrilateralElement, mfem::L2_QuadrilateralElement, mfem::L2Pos_SegmentElement, mfem::L2_SegmentElement, mfem::H1_TetrahedronElement, mfem::H1_TriangleElement, mfem::H1Pos_HexahedronElement, mfem::H1Pos_QuadrilateralElement, mfem::H1Pos_SegmentElement, mfem::H1_HexahedronElement, mfem::H1_QuadrilateralElement, mfem::H1_SegmentElement, mfem::RotTriLinearHexFiniteElement, mfem::RefinedTriLinear3DFiniteElement, mfem::RefinedBiLinear2DFiniteElement, mfem::RefinedLinear3DFiniteElement, mfem::RefinedLinear2DFiniteElement, mfem::RefinedLinear1DFiniteElement, mfem::LagrangeHexFiniteElement, mfem::P0HexFiniteElement, mfem::P0TetFiniteElement, mfem::P1TetNonConfFiniteElement, mfem::Lagrange1DFiniteElement, mfem::P2SegmentFiniteElement, mfem::P1SegmentFiniteElement, mfem::P0SegmentFiniteElement, mfem::CrouzeixRaviartQuadFiniteElement, mfem::CrouzeixRaviartFiniteElement, mfem::TriLinear3DFiniteElement, mfem::Quadratic3DFiniteElement, mfem::Linear3DFiniteElement, mfem::P0QuadFiniteElement, mfem::P0TriangleFiniteElement, mfem::Cubic3DFiniteElement, mfem::Cubic2DFiniteElement, mfem::Cubic1DFiniteElement, mfem::BiCubic2DFiniteElement, mfem::GaussBiQuad2DFiniteElement, mfem::BiQuadPos2DFiniteElement, mfem::BiQuad2DFiniteElement, mfem::GaussQuad2DFiniteElement, mfem::Quad2DFiniteElement, mfem::QuadPos1DFiniteElement, mfem::Quad1DFiniteElement, mfem::P1OnQuadFiniteElement, mfem::GaussBiLinear2DFiniteElement, mfem::GaussLinear2DFiniteElement, mfem::BiLinear2DFiniteElement, mfem::Linear2DFiniteElement, mfem::Linear1DFiniteElement, and mfem::PointFiniteElement.

void mfem::FiniteElement::CalcVShape ( const IntegrationPoint ip,
DenseMatrix shape 
) const
virtual
void mfem::FiniteElement::CalcVShape ( ElementTransformation Trans,
DenseMatrix shape 
) const
virtual
int mfem::FiniteElement::GetDim ( ) const
inline

Returns the space dimension for the finite element.

Definition at line 82 of file fe.hpp.

int mfem::FiniteElement::GetDof ( ) const
inline

Returns the degrees of freedom in the FE space.

Definition at line 88 of file fe.hpp.

void mfem::FiniteElement::GetFaceDofs ( int  face,
int **  dofs,
int *  ndofs 
) const
virtual

Reimplemented in mfem::Linear3DFiniteElement.

Definition at line 58 of file fe.cpp.

int mfem::FiniteElement::GetGeomType ( ) const
inline

Returns the geometry type:

Definition at line 85 of file fe.hpp.

void mfem::FiniteElement::GetLocalInterpolation ( ElementTransformation Trans,
DenseMatrix I 
) const
virtual
int mfem::FiniteElement::GetMapType ( ) const
inline

Definition at line 98 of file fe.hpp.

const IntegrationRule& mfem::FiniteElement::GetNodes ( ) const
inline

Definition at line 113 of file fe.hpp.

int mfem::FiniteElement::GetOrder ( ) const
inline

Returns the order of the finite element.

Definition at line 91 of file fe.hpp.

int mfem::FiniteElement::GetRangeType ( ) const
inline

Definition at line 96 of file fe.hpp.

void mfem::FiniteElement::Project ( Coefficient coeff,
ElementTransformation Trans,
Vector dofs 
) const
virtual

Given a coefficient and a transformation, compute its projection (approximation) in the local finite dimensional space in terms of the degrees of freedom.

Reimplemented in mfem::BiQuadPos2DFiniteElement, mfem::PositiveFiniteElement, and mfem::NodalFiniteElement.

Definition at line 75 of file fe.cpp.

void mfem::FiniteElement::Project ( VectorCoefficient vc,
ElementTransformation Trans,
Vector dofs 
) const
virtual
void mfem::FiniteElement::Project ( const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix I 
) const
virtual
void mfem::FiniteElement::ProjectCurl ( const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix curl 
) const
virtual

Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.

Reimplemented in mfem::ND_TetrahedronElement, mfem::ND_HexahedronElement, mfem::RT_TetrahedronElement, mfem::RT_HexahedronElement, mfem::L2_TriangleElement, and mfem::L2_QuadrilateralElement.

Definition at line 108 of file fe.cpp.

void mfem::FiniteElement::ProjectDelta ( int  vertex,
Vector dofs 
) const
virtual
void mfem::FiniteElement::ProjectDiv ( const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix div 
) const
virtual

Compute the discrete divergence matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.

Reimplemented in mfem::NodalFiniteElement.

Definition at line 116 of file fe.cpp.

void mfem::FiniteElement::ProjectGrad ( const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix grad 
) const
virtual

Compute the discrete gradient matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.

Reimplemented in mfem::ND_SegmentElement, mfem::ND_TriangleElement, mfem::ND_TetrahedronElement, mfem::ND_QuadrilateralElement, mfem::ND_HexahedronElement, mfem::RT_TriangleElement, mfem::RT_QuadrilateralElement, and mfem::NodalFiniteElement.

Definition at line 100 of file fe.cpp.

int mfem::FiniteElement::Space ( ) const
inline

Returns the type of space on each element.

Definition at line 94 of file fe.hpp.

Member Data Documentation

int mfem::FiniteElement::Dim
protected

Definition at line 47 of file fe.hpp.

int mfem::FiniteElement::Dof
protected

Definition at line 47 of file fe.hpp.

int mfem::FiniteElement::FuncSpace
protected

Definition at line 47 of file fe.hpp.

int mfem::FiniteElement::GeomType
protected

Definition at line 47 of file fe.hpp.

int mfem::FiniteElement::MapType
protected

Definition at line 47 of file fe.hpp.

IntegrationRule mfem::FiniteElement::Nodes
protected

Definition at line 48 of file fe.hpp.

int mfem::FiniteElement::Order
protected

Definition at line 47 of file fe.hpp.

int mfem::FiniteElement::RangeType
protected

Definition at line 47 of file fe.hpp.


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