MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Types | Public Member Functions | Static 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 and DerivRangeType. More...
 
enum  { VALUE, INTEGRAL, H_DIV, H_CURL }
 Enumeration for MapType: defines how reference functions are mapped to physical space. More...
 
enum  { NONE, GRAD, DIV, CURL }
 Enumeration for DerivType: defines which derivative method is implemented. More...
 

Public Member Functions

 FiniteElement (int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
 
int GetDim () const
 Returns the reference space dimension for the finite element. More...
 
Geometry::Type GetGeomType () const
 Returns the Geometry::Type of the reference element. More...
 
int GetDof () const
 Returns the number of degrees of freedom in the finite element. More...
 
int GetOrder () const
 Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order. More...
 
bool HasAnisotropicOrders () const
 Returns true if the FiniteElement basis may be using different orders/degrees in different spatial directions. More...
 
const int * GetAnisotropicOrders () const
 Returns an array containing the anisotropic orders/degrees. More...
 
int Space () const
 Returns the type of space on each element. More...
 
int GetRangeType () const
 
int GetDerivRangeType () const
 
int GetMapType () const
 
int GetDerivType () const
 
int GetDerivMapType () const
 
virtual void CalcShape (const IntegrationPoint &ip, Vector &shape) const =0
 Evaluate the values of all shape functions of a scalar finite element in reference space at the given point ip. More...
 
void CalcPhysShape (ElementTransformation &Trans, Vector &shape) const
 Evaluate the values of all shape functions of a scalar finite element in physical space at the point described by Trans. More...
 
virtual void CalcDShape (const IntegrationPoint &ip, DenseMatrix &dshape) const =0
 Evaluate the gradients of all shape functions of a scalar finite element in reference space at the given point ip. More...
 
void CalcPhysDShape (ElementTransformation &Trans, DenseMatrix &dshape) const
 Evaluate the gradients of all shape functions of a scalar finite element in physical space at the point described by Trans. More...
 
const IntegrationRuleGetNodes () const
 
virtual void CalcVShape (const IntegrationPoint &ip, DenseMatrix &shape) const
 Evaluate the values of all shape functions of a vector finite element in reference space at the given point ip. More...
 
virtual void CalcVShape (ElementTransformation &Trans, DenseMatrix &shape) const
 Evaluate the values of all shape functions of a vector finite element in physical space at the point described by Trans. More...
 
void CalcPhysVShape (ElementTransformation &Trans, DenseMatrix &shape) const
 Equivalent to the CalcVShape() method with the same arguments. More...
 
virtual void CalcDivShape (const IntegrationPoint &ip, Vector &divshape) const
 Evaluate the divergence of all shape functions of a vector finite element in reference space at the given point ip. More...
 
void CalcPhysDivShape (ElementTransformation &Trans, Vector &divshape) const
 Evaluate the divergence of all shape functions of a vector finite element in physical space at the point described by Trans. More...
 
virtual void CalcCurlShape (const IntegrationPoint &ip, DenseMatrix &curl_shape) const
 Evaluate the curl of all shape functions of a vector finite element in reference space at the given point ip. More...
 
void CalcPhysCurlShape (ElementTransformation &Trans, DenseMatrix &curl_shape) const
 Evaluate the curl of all shape functions of a vector finite element in physical space at the point described by Trans. More...
 
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
 Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base geometry under the given transformation. More...
 
virtual void GetLocalRestriction (ElementTransformation &Trans, DenseMatrix &R) const
 Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. More...
 
virtual void GetTransferMatrix (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
 Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this finite element. More...
 
virtual void Project (Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
 
virtual void Project (VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
 
virtual void ProjectMatrixCoefficient (MatrixCoefficient &mc, ElementTransformation &T, 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 const DofToQuadGetDofToQuad (const IntegrationRule &ir, DofToQuad::Mode mode) const
 
virtual ~FiniteElement ()
 

Static Public Member Functions

static bool IsClosedType (int b_type)
 
static bool IsOpenType (int b_type)
 
static int VerifyClosed (int b_type)
 
static int VerifyOpen (int b_type)
 
static int VerifyNodal (int b_type)
 

Protected Attributes

int Dim
 Dimension of reference space. More...
 
Geometry::Type GeomType
 Geometry::Type of the reference element. More...
 
int FuncSpace
 
int RangeType
 
int MapType
 
int DerivType
 
int DerivRangeType
 
int DerivMapType
 
int Dof
 Number of degrees of freedom. More...
 
int Order
 Order/degree of the shape functions. More...
 
int Orders [Geometry::MaxDim]
 Anisotropic orders. More...
 
IntegrationRule Nodes
 
DenseMatrix vshape
 
Array< DofToQuad * > dof2quad_array
 Container for all DofToQuad objects created by the FiniteElement. More...
 

Detailed Description

Abstract class for Finite Elements.

Definition at line 229 of file fe.hpp.

Member Enumeration Documentation

anonymous enum

Enumeration for DerivType: defines which derivative method is implemented.

Each FiniteElement class implements only one type of derivative. The value returned by GetDerivType() indicates which derivative method is implemented.

Enumerator
NONE 

No derivatives implemented.

GRAD 

Implements CalcDShape methods.

DIV 

Implements CalcDivShape methods.

CURL 

Implements CalcCurlShape methods.

Definition at line 288 of file fe.hpp.

anonymous enum

Enumeration for RangeType and DerivRangeType.

Enumerator
SCALAR 
VECTOR 

Definition at line 251 of file fe.hpp.

anonymous enum

Enumeration for MapType: defines how reference functions are mapped to physical space.

A reference function, uh(xh), can be mapped to a function, u(x), on a general physical element in following ways:

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 

For scalar fields; preserves point values.

INTEGRAL 

For scalar fields; preserves volume integrals.

H_DIV 

For vector fields; preserves surface integrals of the normal component

H_CURL 

For vector fields; preserves line integrals of the tangential component

Definition at line 273 of file fe.hpp.

Constructor & Destructor Documentation

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

Construct FiniteElement with given

Parameters
DReference space dimension
GGeometry type (of type Geometry::Type)
DoNumber of degrees of freedom in the FiniteElement
OOrder/degree of the FiniteElement
FFunctionSpace type of the FiniteElement

Definition at line 25 of file fe.cpp.

mfem::FiniteElement::~FiniteElement ( )
virtual

Definition at line 214 of file fe.cpp.

Member Function Documentation

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

Evaluate the curl of all shape functions of a vector finite element in reference space at the given point ip.

Each row of the result DenseMatrix curl_shape contains the components of the curl of one vector shape function. The size (Dof x CDim) of curl_shape must be set in advance, where CDim = 3 for Dim = 3 and CDim = 1 for Dim = 2.

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

Definition at line 68 of file fe.cpp.

void mfem::FiniteElement::CalcDivShape ( const IntegrationPoint ip,
Vector divshape 
) const
virtual

Evaluate the divergence of all shape functions of a vector finite element in reference space at the given point ip.

The size (Dof) of the result Vector divshape must be set in advance.

Reimplemented in mfem::RT_TetrahedronElement, mfem::RT_TriangleElement, mfem::RT_HexahedronElement, mfem::RT_QuadrilateralElement, mfem::RT0TetFiniteElement, mfem::RT1HexFiniteElement, mfem::RT0HexFiniteElement, mfem::RT2QuadFiniteElement, mfem::RT2TriangleFiniteElement, mfem::RT1QuadFiniteElement, mfem::RT1TriangleFiniteElement, mfem::RT0QuadFiniteElement, and mfem::RT0TriangleFiniteElement.

Definition at line 54 of file fe.cpp.

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

Evaluate the gradients of all shape functions of a scalar finite element in reference space at the given point ip.

Each row of the result DenseMatrix dshape contains the derivatives of one shape function. The size (Dof x Dim) of dshape must be set in advance.

Implemented in mfem::NURBS3DFiniteElement, mfem::NURBS2DFiniteElement, mfem::NURBS1DFiniteElement, mfem::L2Pos_WedgeElement, mfem::L2_WedgeElement, 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::H1Pos_WedgeElement, mfem::H1_WedgeElement, mfem::H1Pos_TetrahedronElement, mfem::H1Pos_TriangleElement, 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::H1_TetrahedronElement, mfem::H1_TriangleElement, mfem::Cubic2DFiniteElement, mfem::BiCubic2DFiniteElement, mfem::Quad2DFiniteElement, and mfem::BiLinear2DFiniteElement.

Definition at line 105 of file fe.cpp.

void mfem::FiniteElement::CalcPhysCurlShape ( ElementTransformation Trans,
DenseMatrix curl_shape 
) const

Evaluate the curl of all shape functions of a vector finite element in physical space at the point described by Trans.

Each row of the result DenseMatrix curl_shape contains the components of the curl of one vector shape function. The size (Dof x CDim) of curl_shape must be set in advance, where CDim = 3 for Dim = 3 and CDim = 1 for Dim = 2.

Definition at line 75 of file fe.cpp.

void mfem::FiniteElement::CalcPhysDivShape ( ElementTransformation Trans,
Vector divshape 
) const

Evaluate the divergence of all shape functions of a vector finite element in physical space at the point described by Trans.

The size (Dof) of the result Vector divshape must be set in advance.

Definition at line 61 of file fe.cpp.

void mfem::FiniteElement::CalcPhysDShape ( ElementTransformation Trans,
DenseMatrix dshape 
) const

Evaluate the gradients of all shape functions of a scalar finite element in physical space at the point described by Trans.

Each row of the result DenseMatrix dshape contains the derivatives of one shape function. The size (Dof x SDim) of dshape must be set in advance, where SDim >= Dim is the physical space dimension as described by Trans.

Definition at line 195 of file fe.cpp.

void mfem::FiniteElement::CalcPhysShape ( ElementTransformation Trans,
Vector shape 
) const

Evaluate the values of all shape functions of a scalar finite element in physical space at the point described by Trans.

The size (Dof) of the result Vector shape must be set in advance.

Definition at line 185 of file fe.cpp.

void mfem::FiniteElement::CalcPhysVShape ( ElementTransformation Trans,
DenseMatrix shape 
) const
inline

Equivalent to the CalcVShape() method with the same arguments.

Definition at line 386 of file fe.hpp.

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

Evaluate the values of all shape functions of a scalar finite element in reference space at the given point ip.

The size (Dof) of the result Vector shape must be set in advance.

Implemented in mfem::NURBS3DFiniteElement, mfem::NURBS2DFiniteElement, mfem::NURBS1DFiniteElement, mfem::ND_SegmentElement, mfem::L2Pos_WedgeElement, mfem::L2_WedgeElement, 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::H1Pos_WedgeElement, mfem::H1_WedgeElement, mfem::H1Pos_TetrahedronElement, mfem::H1Pos_TriangleElement, 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

Evaluate the values of all shape functions of a vector finite element in physical space at the point described by Trans.

Each row of the result DenseMatrix shape contains the components of one vector shape function. The size (Dof x SDim) of shape must be set in advance, where SDim >= Dim is the physical space dimension as described by Trans.

Reimplemented in mfem::ND_SegmentElement, mfem::ND_TriangleElement, mfem::ND_TetrahedronElement, mfem::ND_QuadrilateralElement, mfem::ND_HexahedronElement, mfem::RT_TetrahedronElement, mfem::RT_TriangleElement, mfem::RT_HexahedronElement, mfem::RT_QuadrilateralElement, mfem::RT0TetFiniteElement, mfem::RT1HexFiniteElement, mfem::RT0HexFiniteElement, mfem::Nedelec1TetFiniteElement, mfem::Nedelec1HexFiniteElement, mfem::RT2QuadFiniteElement, mfem::RT2TriangleFiniteElement, mfem::RT1QuadFiniteElement, mfem::RT1TriangleFiniteElement, mfem::RT0QuadFiniteElement, and mfem::RT0TriangleFiniteElement.

Definition at line 47 of file fe.cpp.

const int* mfem::FiniteElement::GetAnisotropicOrders ( ) const
inline

Returns an array containing the anisotropic orders/degrees.

Definition at line 322 of file fe.hpp.

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

Definition at line 335 of file fe.hpp.

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

Definition at line 329 of file fe.hpp.

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

Definition at line 333 of file fe.hpp.

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

Returns the reference space dimension for the finite element.

Definition at line 305 of file fe.hpp.

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

Returns the number of degrees of freedom in the finite element.

Definition at line 311 of file fe.hpp.

const DofToQuad & mfem::FiniteElement::GetDofToQuad ( const IntegrationRule ir,
DofToQuad::Mode  mode 
) const
virtual

Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mode. See the documentation for DofToQuad for more details.

Reimplemented in mfem::PositiveTensorFiniteElement, mfem::NodalTensorFiniteElement, and mfem::ScalarFiniteElement.

Definition at line 206 of file fe.cpp.

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

Reimplemented in mfem::Linear3DFiniteElement.

Definition at line 100 of file fe.cpp.

Geometry::Type mfem::FiniteElement::GetGeomType ( ) const
inline

Returns the Geometry::Type of the reference element.

Definition at line 308 of file fe.hpp.

void mfem::FiniteElement::GetLocalInterpolation ( ElementTransformation Trans,
DenseMatrix I 
) const
virtual
void mfem::FiniteElement::GetLocalRestriction ( ElementTransformation Trans,
DenseMatrix R 
) const
virtual

Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs.

The fine element is the image of the base geometry under the given transformation, Trans.

The assumption in this method is that a subset of the coarse dofs can be expressed only in terms of the dofs of the given fine element.

Rows in R corresponding to coarse dofs that cannot be expressed in terms of the fine dofs will be marked as invalid by setting the first entry (column 0) in the row to infinity().

This method assumes that the dimensions of R are set before it is called.

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

Definition at line 117 of file fe.cpp.

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

Definition at line 331 of file fe.hpp.

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

Definition at line 364 of file fe.hpp.

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

Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order.

Definition at line 315 of file fe.hpp.

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

Definition at line 327 of file fe.hpp.

void mfem::FiniteElement::GetTransferMatrix ( const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix I 
) const
virtual

Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this finite element.

Trans represents the mapping from the reference element of this element into a subset of the reference space of the element fe, thus allowing the "coarse" FiniteElement to be different from the "fine" FiniteElement as when h-refinement is combined with p-refinement or p-derefinement. It is assumed that both finite elements use the same MapType.

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

Definition at line 123 of file fe.cpp.

bool mfem::FiniteElement::HasAnisotropicOrders ( ) const
inline

Returns true if the FiniteElement basis may be using different orders/degrees in different spatial directions.

Definition at line 319 of file fe.hpp.

static bool mfem::FiniteElement::IsClosedType ( int  b_type)
inlinestatic

Definition at line 521 of file fe.hpp.

static bool mfem::FiniteElement::IsOpenType ( int  b_type)
inlinestatic

Definition at line 528 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 130 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_TriangleElement, mfem::RT_HexahedronElement, mfem::RT_QuadrilateralElement, mfem::L2_TriangleElement, and mfem::L2_QuadrilateralElement.

Definition at line 169 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 177 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 161 of file fe.cpp.

void mfem::FiniteElement::ProjectMatrixCoefficient ( MatrixCoefficient mc,
ElementTransformation T,
Vector dofs 
) const
virtual

Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local finite dimensional space in terms of the degrees of freedom. For VectorFiniteElements, the rows of the coefficient are projected in the vector space.

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

Definition at line 142 of file fe.cpp.

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

Returns the type of space on each element.

Definition at line 325 of file fe.hpp.

static int mfem::FiniteElement::VerifyClosed ( int  b_type)
inlinestatic

Definition at line 535 of file fe.hpp.

static int mfem::FiniteElement::VerifyNodal ( int  b_type)
inlinestatic

Definition at line 546 of file fe.hpp.

static int mfem::FiniteElement::VerifyOpen ( int  b_type)
inlinestatic

Definition at line 541 of file fe.hpp.

Member Data Documentation

int mfem::FiniteElement::DerivMapType
protected

Definition at line 234 of file fe.hpp.

int mfem::FiniteElement::DerivRangeType
protected

Definition at line 234 of file fe.hpp.

int mfem::FiniteElement::DerivType
protected

Definition at line 234 of file fe.hpp.

int mfem::FiniteElement::Dim
protected

Dimension of reference space.

Definition at line 232 of file fe.hpp.

int mfem::FiniteElement::Dof
mutableprotected

Number of degrees of freedom.

Definition at line 237 of file fe.hpp.

Array<DofToQuad*> mfem::FiniteElement::dof2quad_array
mutableprotected

Container for all DofToQuad objects created by the FiniteElement.

Multiple DofToQuad objects may be needed when different quadrature rules or different DofToQuad::Mode are used.

Definition at line 247 of file fe.hpp.

int mfem::FiniteElement::FuncSpace
protected

Definition at line 234 of file fe.hpp.

Geometry::Type mfem::FiniteElement::GeomType
protected

Geometry::Type of the reference element.

Definition at line 233 of file fe.hpp.

int mfem::FiniteElement::MapType
protected

Definition at line 234 of file fe.hpp.

IntegrationRule mfem::FiniteElement::Nodes
protected

Definition at line 240 of file fe.hpp.

int mfem::FiniteElement::Order
mutableprotected

Order/degree of the shape functions.

Definition at line 237 of file fe.hpp.

int mfem::FiniteElement::Orders[Geometry::MaxDim]
mutableprotected

Anisotropic orders.

Definition at line 239 of file fe.hpp.

int mfem::FiniteElement::RangeType
protected

Definition at line 234 of file fe.hpp.

DenseMatrix mfem::FiniteElement::vshape
mutableprotected

Definition at line 242 of file fe.hpp.


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