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

Collection of finite elements from the same family in multiple dimensions. This class is used to match the degrees of freedom of a FiniteElementSpace between elements, and to provide the finite element restriction from an element to its boundary. More...

#include <fe_coll.hpp>

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

Public Types

enum  { CONTINUOUS, TANGENTIAL, NORMAL, DISCONTINUOUS }
 Enumeration for ContType: defines the continuity of the field across element interfaces. More...
 

Public Member Functions

virtual const FiniteElementFiniteElementForGeometry (Geometry::Type GeomType) const =0
 
virtual const FiniteElementFiniteElementForDim (int dim) const
 Returns the first non-NULL FiniteElement for the given dimension. More...
 
virtual int DofForGeometry (Geometry::Type GeomType) const =0
 
virtual StatelessDofTransformationDofTransformationForGeometry (Geometry::Type GeomType) const
 Returns a DoF transformation object compatible with this basis and geometry type. More...
 
virtual const int * DofOrderForOrientation (Geometry::Type GeomType, int Or) const =0
 Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]. More...
 
virtual const char * Name () const
 
virtual int GetContType () const =0
 
int HasFaceDofs (Geometry::Type geom, int p) const
 
virtual const FiniteElementTraceFiniteElementForGeometry (Geometry::Type GeomType) const
 
virtual FiniteElementCollectionGetTraceCollection () const
 
virtual ~FiniteElementCollection ()
 
void SubDofOrder (Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
 Get the local dofs for a given sub-manifold. More...
 
const FiniteElementGetFE (Geometry::Type geom, int p) const
 Variable order version of FiniteElementForGeometry(). More...
 
int GetNumDof (Geometry::Type geom, int p) const
 Variable order version of DofForGeometry(). More...
 
const int * GetDofOrdering (Geometry::Type geom, int p, int ori) const
 Variable order version of DofOrderForOrientation(). More...
 
int GetOrder () const
 Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned by FiniteElement::GetOrder() of the highest-dimensional FiniteElements defined by the collection. More...
 
virtual FiniteElementCollectionClone (int p) const
 Instantiate a new collection of the same type with a different order. More...
 
virtual int GetRangeType (int dim) const
 
virtual int GetDerivRangeType (int dim) const
 
virtual int GetMapType (int dim) const
 
virtual int GetDerivType (int dim) const
 
virtual int GetDerivMapType (int dim) const
 

Static Public Member Functions

static FiniteElementCollectionNew (const char *name)
 Factory method: return a newly allocated FiniteElementCollection according to the given name. More...
 

Protected Types

enum  ErrorMode { RETURN_NULL, RAISE_MFEM_ERROR }
 How to treat errors in FiniteElementForGeometry() calls. More...
 

Protected Member Functions

 FiniteElementCollection ()
 
 FiniteElementCollection (int p)
 
void InitVarOrder (int p) const
 

Static Protected Member Functions

template<Geometry::Type geom>
static void GetNVE (int &nv, int &ne)
 
template<Geometry::Type geom, typename v_t >
static void GetEdge (int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
 
template<Geometry::Type geom, Geometry::Type f_geom, typename v_t , typename e_t , typename eo_t >
static void GetFace (int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, Geometry::Type &fg, int &fo, const int face_info)
 

Protected Attributes

const int base_p
 Order as returned by GetOrder(). More...
 
Array< FiniteElementCollection * > var_orders
 
ErrorMode error_mode = RAISE_MFEM_ERROR
 How to treat errors in FiniteElementForGeometry() calls. More...
 

Detailed Description

Collection of finite elements from the same family in multiple dimensions. This class is used to match the degrees of freedom of a FiniteElementSpace between elements, and to provide the finite element restriction from an element to its boundary.

Definition at line 26 of file fe_coll.hpp.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum

Enumeration for ContType: defines the continuity of the field across element interfaces.

Enumerator
CONTINUOUS 

Field is continuous across element interfaces.

TANGENTIAL 

Tangential components of vector field.

NORMAL 

Normal component of vector field.

DISCONTINUOUS 

Field is discontinuous across element interfaces.

Definition at line 45 of file fe_coll.hpp.

◆ ErrorMode

How to treat errors in FiniteElementForGeometry() calls.

Enumerator
RETURN_NULL 

Return NULL on errors.

RAISE_MFEM_ERROR 

Raise an MFEM error (default in base class). Sub-classes can ignore this and return NULL.

Definition at line 245 of file fe_coll.hpp.

Constructor & Destructor Documentation

◆ ~FiniteElementCollection()

mfem::FiniteElementCollection::~FiniteElementCollection ( )
virtual

Definition at line 378 of file fe_coll.cpp.

◆ FiniteElementCollection() [1/2]

mfem::FiniteElementCollection::FiniteElementCollection ( )
inlineprotected

Definition at line 237 of file fe_coll.hpp.

◆ FiniteElementCollection() [2/2]

mfem::FiniteElementCollection::FiniteElementCollection ( int  p)
inlineprotected

Definition at line 238 of file fe_coll.hpp.

Member Function Documentation

◆ Clone()

FiniteElementCollection * mfem::FiniteElementCollection::Clone ( int  p) const
virtual

Instantiate a new collection of the same type with a different order.

Generally, the order parameter p is NOT the same as the parameter p used by some of the constructors of derived classes. Instead, this p represents the order of the new FE collection as it will be returned by its GetOrder() method.

Reimplemented in mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.

Definition at line 361 of file fe_coll.cpp.

◆ DofForGeometry()

virtual int mfem::FiniteElementCollection::DofForGeometry ( Geometry::Type  GeomType) const
pure virtual

◆ DofOrderForOrientation()

virtual const int* mfem::FiniteElementCollection::DofOrderForOrientation ( Geometry::Type  GeomType,
int  Or 
) const
pure virtual

◆ DofTransformationForGeometry()

virtual StatelessDofTransformation* mfem::FiniteElementCollection::DofTransformationForGeometry ( Geometry::Type  GeomType) const
inlinevirtual

Returns a DoF transformation object compatible with this basis and geometry type.

Reimplemented in mfem::ND_FECollection.

Definition at line 68 of file fe_coll.hpp.

◆ FiniteElementForDim()

const FiniteElement * mfem::FiniteElementCollection::FiniteElementForDim ( int  dim) const
virtual

Returns the first non-NULL FiniteElement for the given dimension.

Note
Repeatedly calls FiniteElementForGeometry in the order defined in the Geometry::Type enumeration.

Definition at line 26 of file fe_coll.cpp.

◆ FiniteElementForGeometry()

virtual const FiniteElement* mfem::FiniteElementCollection::FiniteElementForGeometry ( Geometry::Type  GeomType) const
pure virtual

◆ GetContType()

virtual int mfem::FiniteElementCollection::GetContType ( ) const
pure virtual

◆ GetDerivMapType()

int mfem::FiniteElementCollection::GetDerivMapType ( int  dim) const
virtual

Definition at line 80 of file fe_coll.cpp.

◆ GetDerivRangeType()

int mfem::FiniteElementCollection::GetDerivRangeType ( int  dim) const
virtual

Definition at line 50 of file fe_coll.cpp.

◆ GetDerivType()

int mfem::FiniteElementCollection::GetDerivType ( int  dim) const
virtual

Definition at line 70 of file fe_coll.cpp.

◆ GetDofOrdering()

const int* mfem::FiniteElementCollection::GetDofOrdering ( Geometry::Type  geom,
int  p,
int  ori 
) const
inline

Variable order version of DofOrderForOrientation().

The order parameter p represents the order of the highest-dimensional FiniteElements the fixed-order collection we want to query. In general, this order is different from the order of the element corresponding to geom in that fixed-order collection.

Definition at line 215 of file fe_coll.hpp.

◆ GetEdge()

template<Geometry::Type geom, typename v_t >
void mfem::FiniteElementCollection::GetEdge ( int &  nv,
v_t &  v,
int &  ne,
int &  e,
int &  eo,
const int  edge_info 
)
inlinestaticprotected

Definition at line 397 of file fe_coll.cpp.

◆ GetFace()

template<Geometry::Type geom, Geometry::Type f_geom, typename v_t , typename e_t , typename eo_t >
void mfem::FiniteElementCollection::GetFace ( int &  nv,
v_t &  v,
int &  ne,
e_t &  e,
eo_t &  eo,
int &  nf,
int &  f,
Geometry::Type fg,
int &  fo,
const int  face_info 
)
inlinestaticprotected

Definition at line 417 of file fe_coll.cpp.

◆ GetFE()

const FiniteElement* mfem::FiniteElementCollection::GetFE ( Geometry::Type  geom,
int  p 
) const
inline

Variable order version of FiniteElementForGeometry().

The order parameter p represents the order of the highest-dimensional FiniteElements the fixed-order collection we want to query. In general, this order is different from the order of the returned FiniteElement.

Definition at line 191 of file fe_coll.hpp.

◆ GetMapType()

int mfem::FiniteElementCollection::GetMapType ( int  dim) const
virtual

Definition at line 60 of file fe_coll.cpp.

◆ GetNumDof()

int mfem::FiniteElementCollection::GetNumDof ( Geometry::Type  geom,
int  p 
) const
inline

Variable order version of DofForGeometry().

The order parameter p represents the order of the highest-dimensional FiniteElements the fixed-order collection we want to query. In general, this order is different from the order of the element corresponding to geom in that fixed-order collection.

Definition at line 203 of file fe_coll.hpp.

◆ GetNVE()

template<Geometry::Type geom>
void mfem::FiniteElementCollection::GetNVE ( int &  nv,
int &  ne 
)
inlinestaticprotected

Definition at line 387 of file fe_coll.cpp.

◆ GetOrder()

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

Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned by FiniteElement::GetOrder() of the highest-dimensional FiniteElements defined by the collection.

Definition at line 225 of file fe_coll.hpp.

◆ GetRangeType()

int mfem::FiniteElementCollection::GetRangeType ( int  dim) const
virtual
Note
The following methods provide the same information as the corresponding methods of the FiniteElement base class.

Definition at line 40 of file fe_coll.cpp.

◆ GetTraceCollection()

FiniteElementCollection * mfem::FiniteElementCollection::GetTraceCollection ( ) const
virtual

◆ HasFaceDofs()

int mfem::FiniteElementCollection::HasFaceDofs ( Geometry::Type  geom,
int  p 
) const

Definition at line 90 of file fe_coll.cpp.

◆ InitVarOrder()

void mfem::FiniteElementCollection::InitVarOrder ( int  p) const
protected

Definition at line 369 of file fe_coll.cpp.

◆ Name()

virtual const char* mfem::FiniteElementCollection::Name ( ) const
inlinevirtual

◆ New()

FiniteElementCollection * mfem::FiniteElementCollection::New ( const char *  name)
static

Factory method: return a newly allocated FiniteElementCollection according to the given name.

FEC Name Space Order BasisType FiniteElement::M
H1_[DIM]_[ORDER] H1 * 1 VALUE H1 nodal elements
H1@[BTYPE]_[DIM]_[ORDER] H1 * * VALUE H1 nodal elements
H1Pos_[DIM]_[ORDER] H1 * 1 VALUE H1 nodal elements
H1Pos_Trace_[DIM]_[ORDER] H^{1/2} * 2 VALUE H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices)
H1_Trace_[DIM]_[ORDER] H^{1/2} * 1 VALUE H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices)
H1_Trace@[BTYPE]_[DIM]_[ORDER] H^{1/2} * 1 VALUE H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices)
ND_[DIM]_[ORDER] H(curl) * 1 / 0 H_CURL Nedelec vector elements
ND@[CBTYPE][OBTYPE]_[DIM]_[ORDER] H(curl) * * / * H_CURL Nedelec vector elements
ND_Trace_[DIM]_[ORDER] H^{1/2} * 1 / 0 H_CURL H^{1/2}-conforming trace elements for H(curl) defined on the interface between mesh elements (faces)
ND_Trace@[CBTYPE][OBTYPE]_[DIM]_[ORDER] H^{1/2} * 1 / 0 H_CURL H^{1/2}-conforming trace elements for H(curl) defined on the interface between mesh elements (faces)
RT_[DIM]_[ORDER] H(div) * 1 / 0 H_DIV Raviart-Thomas vector elements
RT@[CBTYPE][OBTYPE]_[DIM]_[ORDER] H(div) * * / * H_DIV Raviart-Thomas vector elements
RT_Trace_[DIM]_[ORDER] H^{1/2} * 1 / 0 INTEGRAL H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces)
RT_ValTrace_[DIM]_[ORDER] H^{1/2} * 1 / 0 VALUE H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces)
RT_Trace@[BTYPE]_[DIM]_[ORDER] H^{1/2} * 1 / 0 INTEGRAL H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces)
RT_ValTrace@[BTYPE]_[DIM]_[ORDER] H^{1/2} * 1 / 0 VALUE H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces)
L2_[DIM]_[ORDER] L2 * 0 VALUE Discontinuous L2 elements
L2_T[BTYPE]_[DIM]_[ORDER] L2 * 0 VALUE Discontinuous L2 elements
L2Int_[DIM]_[ORDER] L2 * 0 INTEGRAL Discontinuous L2 elements
L2Int_T[BTYPE]_[DIM]_[ORDER] L2 * 0 INTEGRAL Discontinuous L2 elements
DG_Iface_[DIM]_[ORDER] - * 0 VALUE Discontinuous elements on the interface between mesh elements (faces)
DG_Iface@[BTYPE]_[DIM]_[ORDER] - * 0 VALUE Discontinuous elements on the interface between mesh elements (faces)
DG_IntIface_[DIM]_[ORDER] - * 0 INTEGRAL Discontinuous elements on the interface between mesh elements (faces)
DG_IntIface@[BTYPE]_[DIM]_[ORDER] - * 0 INTEGRAL Discontinuous elements on the interface between mesh elements (faces)
NURBS[ORDER] - * - VALUE Non-Uniform Rational B-Splines (NURBS) elements
LinearNonConf3D - 1 1 VALUE Piecewise-linear nonconforming finite elements in 3D
CrouzeixRaviart - - - - Crouzeix-Raviart nonconforming elements in 2D
Local_[FENAME] - - - - Special collection that builds a local version out of the FENAME collection
- - - - - -
Linear H1 1 1 VALUE Left in for backward compatibility, consider using H1_
Quadratic H1 2 1 VALUE Left in for backward compatibility, consider using H1_
QuadraticPos H1 2 2 VALUE Left in for backward compatibility, consider using H1_
Cubic H1 2 1 VALUE Left in for backward compatibility, consider using H1_
Const2D L2 0 1 VALUE Left in for backward compatibility, consider using L2_
Const3D L2 0 1 VALUE Left in for backward compatibility, consider using L2_
LinearDiscont2D L2 1 1 VALUE Left in for backward compatibility, consider using L2_
GaussLinearDiscont2D L2 1 0 VALUE Left in for backward compatibility, consider using L2_
P1OnQuad H1 1 1 VALUE Linear P1 element with 3 nodes on a square
QuadraticDiscont2D L2 2 1 VALUE Left in for backward compatibility, consider using L2_
QuadraticPosDiscont2D L2 2 2 VALUE Left in for backward compatibility, consider using L2_
GaussQuadraticDiscont2D L2 2 0 VALUE Left in for backward compatibility, consider using L2_
CubicDiscont2D L2 3 1 VALUE Left in for backward compatibility, consider using L2_
LinearDiscont3D L2 1 1 VALUE Left in for backward compatibility, consider using L2_
QuadraticDiscont3D L2 2 1 VALUE Left in for backward compatibility, consider using L2_
ND1_3D H(Curl) 1 1 / 0 H_CURL Left in for backward compatibility, consider using ND_
RT0_2D H(Div) 1 1 / 0 H_DIV Left in for backward compatibility, consider using RT_
RT1_2D H(Div) 2 1 / 0 H_DIV Left in for backward compatibility, consider using RT_
RT2_2D H(Div) 3 1 / 0 H_DIV Left in for backward compatibility, consider using RT_
RT0_3D H(Div) 1 1 / 0 H_DIV Left in for backward compatibility, consider using RT_
RT1_3D H(Div) 2 1 / 0 H_DIV Left in for backward compatibility, consider using RT_
Tag Description
[DIM] Dimension of the elements (1D, 2D, 3D)
[ORDER] Approximation order of the elements (P0, P1, P2, ...)
[BTYPE] BasisType of the element (0-GaussLegendre, 1 - GaussLobatto, 2-Bernstein, 3-OpenUniform, 4-CloseUniform, 5-OpenHalfUniform)
[OBTYPE] Open BasisType of the element for elements which have both types
[CBTYPE] Closed BasisType of the element for elements which have both types

[FENAME] Is a special case for the Local FEC which generates a local version of a given FEC. It is selected from one of (BiCubic2DFiniteElement, Quad_Q3, Nedelec1HexFiniteElement, Hex_ND1, H1_[DIM]_[ORDER],H1Pos_[DIM]_[ORDER], L2_[DIM]_[ORDER] )

Definition at line 116 of file fe_coll.cpp.

◆ SubDofOrder()

void mfem::FiniteElementCollection::SubDofOrder ( Geometry::Type  Geom,
int  SDim,
int  Info,
Array< int > &  dofs 
) const

Get the local dofs for a given sub-manifold.

Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex, 1D

  • edge, 2D - face) including those on its boundary. The local index of the sub-manifold (inside Geom) and its orientation are given by the parameter Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed that 0 <= SDim <= Dim(Geom).

Definition at line 462 of file fe_coll.cpp.

◆ TraceFiniteElementForGeometry()

virtual const FiniteElement* mfem::FiniteElementCollection::TraceFiniteElementForGeometry ( Geometry::Type  GeomType) const
inlinevirtual

Reimplemented in mfem::L2_FECollection.

Definition at line 97 of file fe_coll.hpp.

Member Data Documentation

◆ base_p

const int mfem::FiniteElementCollection::base_p
protected

Order as returned by GetOrder().

Definition at line 235 of file fe_coll.hpp.

◆ error_mode

ErrorMode mfem::FiniteElementCollection::error_mode = RAISE_MFEM_ERROR
mutableprotected

How to treat errors in FiniteElementForGeometry() calls.

The typical error in derived classes is that no FiniteElement is defined for the given Geometry, or the input is not a valid Geometry.

Definition at line 255 of file fe_coll.hpp.

◆ var_orders

Array<FiniteElementCollection*> mfem::FiniteElementCollection::var_orders
mutableprotected

Definition at line 242 of file fe_coll.hpp.


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