MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
mfem::VectorFiniteElement Class Reference

Intermediate class for finite elements whose basis functions return vector values. More...

#include <fe_base.hpp>

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

Public Member Functions

 VectorFiniteElement (int D, Geometry::Type G, int Do, int O, int M, int F=FunctionSpace::Pk)
 
- Public Member Functions inherited from mfem::FiniteElement
 FiniteElement (int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
 Construct FiniteElement with given. More...
 
int GetDim () const
 Returns the reference space dimension for the finite element. More...
 
int GetVDim () const
 Returns the vector dimension for vector-valued finite elements. More...
 
int GetCurlDim () const
 Returns the dimension of the curl for vector-valued finite elements. 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 FunctionSpace on the element. More...
 
int GetRangeType () const
 Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}. More...
 
int GetDerivRangeType () const
 Returns the FiniteElement::RangeType of the element derivative, either SCALAR or VECTOR. More...
 
int GetMapType () const
 Returns the FiniteElement::MapType of the element describing how reference functions are mapped to physical space, one of {VALUE, INTEGRAL H_DIV, H_CURL}. More...
 
int GetDerivType () const
 Returns the FiniteElement::DerivType of the element describing the spatial derivative method implemented, one of {NONE, GRAD, DIV, CURL}. More...
 
int GetDerivMapType () const
 Returns the FiniteElement::DerivType of the element describing how reference function derivatives are mapped to physical space, one of {VALUE, INTEGRAL, H_DIV, H_CURL}. 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...
 
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
 Get a const reference to the nodes of the element. More...
 
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...
 
virtual 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
 Get the dofs associated with the given face. *dofs is set to an internal array of the local dofc on the face, while *ndofs is set to the number of dofs on that face. More...
 
virtual void CalcHessian (const IntegrationPoint &ip, DenseMatrix &Hessian) const
 Evaluate the Hessians of all shape functions of a scalar finite element in reference space at the given point ip. More...
 
virtual void CalcPhysHessian (ElementTransformation &Trans, DenseMatrix &Hessian) const
 Evaluate the Hessian of all shape functions of a scalar finite element in reference space at the given point ip. More...
 
virtual void CalcPhysLaplacian (ElementTransformation &Trans, Vector &Laplacian) const
 Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the given point ip. More...
 
virtual void CalcPhysLinLaplacian (ElementTransformation &Trans, Vector &Laplacian) 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
 Given a coefficient and a transformation, compute its projection (approximation) in the local finite dimensional space in terms of the degrees of freedom. More...
 
virtual void Project (VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
 Given a vector coefficient and a transformation, compute its projection (approximation) in the local finite dimensional space in terms of the degrees of freedom. (VectorFiniteElements) More...
 
virtual void ProjectFromNodes (Vector &vc, ElementTransformation &Trans, Vector &dofs) const
 Given a vector of values at the finite element nodes and a transformation, compute its projection (approximation) in the local finite dimensional space in terms of the degrees of freedom. Valid for VectorFiniteElements. More...
 
virtual void ProjectMatrixCoefficient (MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
 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. More...
 
virtual void ProjectDelta (int vertex, Vector &dofs) const
 Project a delta function centered on the given vertex in the local finite dimensional space represented by the dofs. More...
 
virtual void Project (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
 Compute the embedding/projection matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the projection depends on it. More...
 
virtual void ProjectGrad (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
 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. More...
 
virtual void ProjectCurl (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
 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. More...
 
virtual void ProjectDiv (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
 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. More...
 
virtual const DofToQuadGetDofToQuad (const IntegrationRule &ir, DofToQuad::Mode mode) const
 Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mode. More...
 
virtual void GetFaceMap (const int face_id, Array< int > &face_map) const
 Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local face face_id. More...
 
virtual StatelessDofTransformationGetDofTransformation () const
 Return a DoF transformation object for this particular type of basis. More...
 
virtual ~FiniteElement ()
 Deconstruct the FiniteElement. More...
 

Protected Member Functions

void SetDerivMembers ()
 
void CalcVShape_RT (ElementTransformation &Trans, DenseMatrix &shape) const
 
void CalcVShape_ND (ElementTransformation &Trans, DenseMatrix &shape) const
 
void Project_RT (const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
 Project a vector coefficient onto the RT basis functions. More...
 
void Project_RT (const double *nk, const Array< int > &d2n, Vector &vc, ElementTransformation &Trans, Vector &dofs) const
 Projects the vector of values given at FE nodes to RT space. More...
 
void ProjectMatrixCoefficient_RT (const double *nk, const Array< int > &d2n, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
 Project the rows of the matrix coefficient in an RT space. More...
 
void Project_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
 Project vector-valued basis functions onto the RT basis functions. More...
 
void ProjectGrad_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
 
void ProjectCurl_ND (const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
 
void ProjectCurl_RT (const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
 
void Project_ND (const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
 Project a vector coefficient onto the ND basis functions. More...
 
void Project_ND (const double *tk, const Array< int > &d2t, Vector &vc, ElementTransformation &Trans, Vector &dofs) const
 Projects the vector of values given at FE nodes to ND space. More...
 
void ProjectMatrixCoefficient_ND (const double *tk, const Array< int > &d2t, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
 Project the rows of the matrix coefficient in an ND space. More...
 
void Project_ND (const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
 Project vector-valued basis functions onto the ND basis functions. More...
 
void ProjectGrad_ND (const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
 
void LocalL2Projection_RT (const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
 
void LocalInterpolation_RT (const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
 
void LocalL2Projection_ND (const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
 
void LocalInterpolation_ND (const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
 
void LocalRestriction_RT (const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
 
void LocalRestriction_ND (const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
 

Static Protected Member Functions

static const VectorFiniteElementCheckVectorFE (const FiniteElement &fe)
 

Protected Attributes

bool is_nodal
 
DenseMatrix JtJ
 
DenseMatrix curlshape
 
DenseMatrix curlshape_J
 
- Protected Attributes inherited from mfem::FiniteElement
int dim
 Dimension of reference space. More...
 
int vdim
 Vector dimension of vector-valued basis functions. More...
 
int cdim
 Dimension of curl for vector-valued basis functions. More...
 
Geometry::Type geom_type
 Geometry::Type of the reference element. More...
 
int func_space
 
int range_type
 
int map_type
 
int deriv_type
 
int deriv_range_type
 
int deriv_map_type
 
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...
 

Additional Inherited Members

- Public Types inherited from mfem::FiniteElement
enum  RangeType { UNKNOWN_RANGE_TYPE = -1, SCALAR, VECTOR }
 Enumeration for range_type and deriv_range_type. More...
 
enum  MapType {
  UNKNOWN_MAP_TYPE = -1, VALUE, INTEGRAL, H_DIV,
  H_CURL
}
 Enumeration for MapType: defines how reference functions are mapped to physical space. More...
 
enum  DerivType { NONE, GRAD, DIV, CURL }
 Enumeration for DerivType: defines which derivative method is implemented. More...
 
- Static Public Member Functions inherited from mfem::FiniteElement
static bool IsClosedType (int b_type)
 Return true if the BasisType of b_type is closed (has Quadrature1D points on the boundary). More...
 
static bool IsOpenType (int b_type)
 Return true if the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary). More...
 
static int VerifyClosed (int b_type)
 Ensure that the BasisType of b_type is closed (has Quadrature1D points on the boundary). More...
 
static int VerifyOpen (int b_type)
 Ensure that the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary). More...
 
static int VerifyNodal (int b_type)
 Ensure that the BasisType of b_type nodal (satisfies the interpolation property). More...
 

Detailed Description

Intermediate class for finite elements whose basis functions return vector values.

Definition at line 788 of file fe_base.hpp.

Constructor & Destructor Documentation

◆ VectorFiniteElement()

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

Definition at line 891 of file fe_base.cpp.

Member Function Documentation

◆ CalcVShape_ND()

void mfem::VectorFiniteElement::CalcVShape_ND ( ElementTransformation Trans,
DenseMatrix shape 
) const
protected

Definition at line 969 of file fe_base.cpp.

◆ CalcVShape_RT()

void mfem::VectorFiniteElement::CalcVShape_RT ( ElementTransformation Trans,
DenseMatrix shape 
) const
protected

Definition at line 957 of file fe_base.cpp.

◆ CheckVectorFE()

static const VectorFiniteElement& mfem::VectorFiniteElement::CheckVectorFE ( const FiniteElement fe)
inlinestaticprotected

Definition at line 949 of file fe_base.hpp.

◆ LocalInterpolation_ND()

void mfem::VectorFiniteElement::LocalInterpolation_ND ( const VectorFiniteElement cfe,
const double *  tk,
const Array< int > &  d2t,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Definition at line 1487 of file fe_base.cpp.

◆ LocalInterpolation_RT()

void mfem::VectorFiniteElement::LocalInterpolation_RT ( const VectorFiniteElement cfe,
const double *  nk,
const Array< int > &  d2n,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Definition at line 1400 of file fe_base.cpp.

◆ LocalL2Projection_ND()

void mfem::VectorFiniteElement::LocalL2Projection_ND ( const VectorFiniteElement cfe,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Definition at line 1441 of file fe_base.cpp.

◆ LocalL2Projection_RT()

void mfem::VectorFiniteElement::LocalL2Projection_RT ( const VectorFiniteElement cfe,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Definition at line 1353 of file fe_base.cpp.

◆ LocalRestriction_ND()

void mfem::VectorFiniteElement::LocalRestriction_ND ( const double *  tk,
const Array< int > &  d2t,
ElementTransformation Trans,
DenseMatrix R 
) const
protected

Definition at line 1569 of file fe_base.cpp.

◆ LocalRestriction_RT()

void mfem::VectorFiniteElement::LocalRestriction_RT ( const double *  nk,
const Array< int > &  d2n,
ElementTransformation Trans,
DenseMatrix R 
) const
protected

Definition at line 1526 of file fe_base.cpp.

◆ Project_ND() [1/3]

void mfem::VectorFiniteElement::Project_ND ( const double *  tk,
const Array< int > &  d2t,
VectorCoefficient vc,
ElementTransformation Trans,
Vector dofs 
) const
protected

Project a vector coefficient onto the ND basis functions.

Parameters
tkEdge tangent vectors for this element type
d2tOffset into tk for each degree of freedom
vcVector coefficient to be projected
TransTransformation from reference to physical coordinates
dofsExpansion coefficients for the approximation of vc

Definition at line 1204 of file fe_base.cpp.

◆ Project_ND() [2/3]

void mfem::VectorFiniteElement::Project_ND ( const double *  tk,
const Array< int > &  d2t,
Vector vc,
ElementTransformation Trans,
Vector dofs 
) const
protected

Projects the vector of values given at FE nodes to ND space.

Project vector values onto the ND basis functions

Parameters
tkEdge tangent vectors for this element type
d2tOffset into tk for each degree of freedom
vcVector values at each interpolation point
TransTransformation from reference to physical coordinates
dofsExpansion coefficients for the approximation of vc

Definition at line 1221 of file fe_base.cpp.

◆ Project_ND() [3/3]

void mfem::VectorFiniteElement::Project_ND ( const double *  tk,
const Array< int > &  d2t,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Project vector-valued basis functions onto the ND basis functions.

Parameters
tkEdge tangent vectors for this element type
d2tOffset into tk for each degree of freedom
feVector-valued finite element basis
TransTransformation from reference to physical coordinates
IExpansion coefficients for the approximation of each basis function

Note: If the FiniteElement, fe, is scalar-valued the projection will assume that a FiniteElementSpace is being used to define a vector field using the scalar basis functions for each component of the vector field.

Definition at line 1259 of file fe_base.cpp.

◆ Project_RT() [1/3]

void mfem::VectorFiniteElement::Project_RT ( const double *  nk,
const Array< int > &  d2n,
VectorCoefficient vc,
ElementTransformation Trans,
Vector dofs 
) const
protected

Project a vector coefficient onto the RT basis functions.

Parameters
nkFace normal vectors for this element type
d2nOffset into nk for each degree of freedom
vcVector coefficient to be projected
TransTransformation from reference to physical coordinates
dofsExpansion coefficients for the approximation of vc

Definition at line 980 of file fe_base.cpp.

◆ Project_RT() [2/3]

void mfem::VectorFiniteElement::Project_RT ( const double *  nk,
const Array< int > &  d2n,
Vector vc,
ElementTransformation Trans,
Vector dofs 
) const
protected

Projects the vector of values given at FE nodes to RT space.

Project vector values onto the RT basis functions

Parameters
nkFace normal vectors for this element type
d2nOffset into nk for each degree of freedom
vcVector values at each interpolation point
TransTransformation from reference to physical coordinates
dofsExpansion coefficients for the approximation of vc

Definition at line 1000 of file fe_base.cpp.

◆ Project_RT() [3/3]

void mfem::VectorFiniteElement::Project_RT ( const double *  nk,
const Array< int > &  d2n,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix I 
) const
protected

Project vector-valued basis functions onto the RT basis functions.

Parameters
nkFace normal vectors for this element type
d2nOffset into nk for each degree of freedom
feVector-valued finite element basis
TransTransformation from reference to physical coordinates
IExpansion coefficients for the approximation of each basis function

Note: If the FiniteElement, fe, is scalar-valued the projection will assume that a FiniteElementSpace is being used to define a vector field using the scalar basis functions for each component of the vector field.

Definition at line 1045 of file fe_base.cpp.

◆ ProjectCurl_ND()

void mfem::VectorFiniteElement::ProjectCurl_ND ( const double *  tk,
const Array< int > &  d2t,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix curl 
) const
protected

Definition at line 1147 of file fe_base.cpp.

◆ ProjectCurl_RT()

void mfem::VectorFiniteElement::ProjectCurl_RT ( const double *  nk,
const Array< int > &  d2n,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix curl 
) const
protected

Definition at line 1185 of file fe_base.cpp.

◆ ProjectGrad_ND()

void mfem::VectorFiniteElement::ProjectGrad_ND ( const double *  tk,
const Array< int > &  d2t,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix grad 
) const
protected

Definition at line 1332 of file fe_base.cpp.

◆ ProjectGrad_RT()

void mfem::VectorFiniteElement::ProjectGrad_RT ( const double *  nk,
const Array< int > &  d2n,
const FiniteElement fe,
ElementTransformation Trans,
DenseMatrix grad 
) const
protected

Definition at line 1120 of file fe_base.cpp.

◆ ProjectMatrixCoefficient_ND()

void mfem::VectorFiniteElement::ProjectMatrixCoefficient_ND ( const double *  tk,
const Array< int > &  d2t,
MatrixCoefficient mc,
ElementTransformation T,
Vector dofs 
) const
protected

Project the rows of the matrix coefficient in an ND space.

Definition at line 1233 of file fe_base.cpp.

◆ ProjectMatrixCoefficient_RT()

void mfem::VectorFiniteElement::ProjectMatrixCoefficient_RT ( const double *  nk,
const Array< int > &  d2n,
MatrixCoefficient mc,
ElementTransformation T,
Vector dofs 
) const
protected

Project the rows of the matrix coefficient in an RT space.

Definition at line 1017 of file fe_base.cpp.

◆ SetDerivMembers()

void mfem::VectorFiniteElement::SetDerivMembers ( )
protected

Definition at line 920 of file fe_base.cpp.

Member Data Documentation

◆ curlshape

DenseMatrix mfem::VectorFiniteElement::curlshape
mutableprotected

Definition at line 804 of file fe_base.hpp.

◆ curlshape_J

DenseMatrix mfem::VectorFiniteElement::curlshape_J
mutableprotected

Definition at line 804 of file fe_base.hpp.

◆ is_nodal

bool mfem::VectorFiniteElement::is_nodal
protected

Definition at line 801 of file fe_base.hpp.

◆ JtJ

DenseMatrix mfem::VectorFiniteElement::JtJ
mutableprotected

Definition at line 803 of file fe_base.hpp.


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