![]() |
MFEM
v4.2.0
Finite element discretization library
|
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>
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 FiniteElement * | FiniteElementForGeometry (Geometry::Type GeomType) const =0 |
| virtual int | DofForGeometry (Geometry::Type GeomType) const =0 |
| 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 GeomType) const |
| virtual const FiniteElement * | TraceFiniteElementForGeometry (Geometry::Type GeomType) const |
| virtual FiniteElementCollection * | GetTraceCollection () 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... | |
Static Public Member Functions | |
| static FiniteElementCollection * | New (const char *name) |
| Factory method: return a newly allocated FiniteElementCollection according to the given name. More... | |
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) |
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.
| anonymous enum |
Enumeration for ContType: defines the continuity of the field across element interfaces.
Definition at line 45 of file fe_coll.hpp.
|
inlinevirtual |
Definition at line 75 of file fe_coll.hpp.
|
pure virtual |
Implemented in mfem::Local_FECollection, mfem::RT1_3DFECollection, mfem::RT0_3DFECollection, mfem::ND1_3DFECollection, mfem::RefinedLinearFECollection, mfem::QuadraticDiscont3DFECollection, mfem::LinearDiscont3DFECollection, mfem::Const3DFECollection, mfem::CubicDiscont2DFECollection, mfem::GaussQuadraticDiscont2DFECollection, mfem::QuadraticPosDiscont2DFECollection, mfem::QuadraticDiscont2DFECollection, mfem::P1OnQuadFECollection, mfem::GaussLinearDiscont2DFECollection, mfem::LinearDiscont2DFECollection, mfem::Const2DFECollection, mfem::RT2_2DFECollection, mfem::RT1_2DFECollection, mfem::RT0_2DFECollection, mfem::LinearNonConf3DFECollection, mfem::CrouzeixRaviartFECollection, mfem::CubicFECollection, mfem::QuadraticPosFECollection, mfem::QuadraticFECollection, mfem::LinearFECollection, mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.
|
pure virtual |
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i].
Implemented in mfem::Local_FECollection, mfem::RT1_3DFECollection, mfem::RT0_3DFECollection, mfem::ND1_3DFECollection, mfem::RefinedLinearFECollection, mfem::QuadraticDiscont3DFECollection, mfem::LinearDiscont3DFECollection, mfem::Const3DFECollection, mfem::CubicDiscont2DFECollection, mfem::GaussQuadraticDiscont2DFECollection, mfem::QuadraticPosDiscont2DFECollection, mfem::QuadraticDiscont2DFECollection, mfem::P1OnQuadFECollection, mfem::GaussLinearDiscont2DFECollection, mfem::LinearDiscont2DFECollection, mfem::Const2DFECollection, mfem::RT2_2DFECollection, mfem::RT1_2DFECollection, mfem::RT0_2DFECollection, mfem::LinearNonConf3DFECollection, mfem::CrouzeixRaviartFECollection, mfem::CubicFECollection, mfem::QuadraticPosFECollection, mfem::QuadraticFECollection, mfem::LinearFECollection, mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.
|
pure virtual |
Implemented in mfem::Local_FECollection, mfem::RT1_3DFECollection, mfem::RT0_3DFECollection, mfem::ND1_3DFECollection, mfem::RefinedLinearFECollection, mfem::QuadraticDiscont3DFECollection, mfem::LinearDiscont3DFECollection, mfem::Const3DFECollection, mfem::CubicDiscont2DFECollection, mfem::GaussQuadraticDiscont2DFECollection, mfem::QuadraticPosDiscont2DFECollection, mfem::QuadraticDiscont2DFECollection, mfem::P1OnQuadFECollection, mfem::GaussLinearDiscont2DFECollection, mfem::LinearDiscont2DFECollection, mfem::Const2DFECollection, mfem::RT2_2DFECollection, mfem::RT1_2DFECollection, mfem::RT0_2DFECollection, mfem::LinearNonConf3DFECollection, mfem::CrouzeixRaviartFECollection, mfem::CubicFECollection, mfem::QuadraticPosFECollection, mfem::QuadraticFECollection, mfem::LinearFECollection, mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.
|
pure virtual |
Implemented in mfem::Local_FECollection, mfem::RT1_3DFECollection, mfem::RT0_3DFECollection, mfem::ND1_3DFECollection, mfem::RefinedLinearFECollection, mfem::QuadraticDiscont3DFECollection, mfem::LinearDiscont3DFECollection, mfem::Const3DFECollection, mfem::CubicDiscont2DFECollection, mfem::GaussQuadraticDiscont2DFECollection, mfem::QuadraticPosDiscont2DFECollection, mfem::QuadraticDiscont2DFECollection, mfem::P1OnQuadFECollection, mfem::GaussLinearDiscont2DFECollection, mfem::LinearDiscont2DFECollection, mfem::Const2DFECollection, mfem::RT2_2DFECollection, mfem::RT1_2DFECollection, mfem::RT0_2DFECollection, mfem::LinearNonConf3DFECollection, mfem::CrouzeixRaviartFECollection, mfem::CubicFECollection, mfem::QuadraticPosFECollection, mfem::QuadraticFECollection, mfem::LinearFECollection, mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.
|
inlinestaticprotected |
Definition at line 303 of file fe_coll.cpp.
|
inlinestaticprotected |
Definition at line 323 of file fe_coll.cpp.
|
inlinestaticprotected |
Definition at line 293 of file fe_coll.cpp.
|
virtual |
Reimplemented in mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, and mfem::H1_FECollection.
Definition at line 41 of file fe_coll.cpp.
| int mfem::FiniteElementCollection::HasFaceDofs | ( | Geometry::Type | GeomType | ) | const |
Definition at line 25 of file fe_coll.cpp.
|
inlinevirtual |
Reimplemented in mfem::Local_FECollection, mfem::RT1_3DFECollection, mfem::RT0_3DFECollection, mfem::ND1_3DFECollection, mfem::RefinedLinearFECollection, mfem::QuadraticDiscont3DFECollection, mfem::LinearDiscont3DFECollection, mfem::Const3DFECollection, mfem::CubicDiscont2DFECollection, mfem::GaussQuadraticDiscont2DFECollection, mfem::QuadraticPosDiscont2DFECollection, mfem::QuadraticDiscont2DFECollection, mfem::P1OnQuadFECollection, mfem::GaussLinearDiscont2DFECollection, mfem::LinearDiscont2DFECollection, mfem::Const2DFECollection, mfem::RT2_2DFECollection, mfem::RT1_2DFECollection, mfem::RT0_2DFECollection, mfem::LinearNonConf3DFECollection, mfem::CrouzeixRaviartFECollection, mfem::CubicFECollection, mfem::QuadraticPosFECollection, mfem::QuadraticFECollection, mfem::LinearFECollection, mfem::NURBSFECollection, mfem::ND_FECollection, mfem::RT_FECollection, mfem::L2_FECollection, and mfem::H1_FECollection.
Definition at line 61 of file fe_coll.hpp.
|
static |
Factory method: return a newly allocated FiniteElementCollection according to the given name.
| FEC Name | Space | Order | BasisType | FiniteElement::MapT | Notes |
|---|---|---|---|---|---|
| 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 | Discontinous L2 elements |
| L2_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinous L2 elements |
| L2Int_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinous L2 elements |
| L2Int_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinous 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 47 of file fe_coll.cpp.
| 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
Definition at line 368 of file fe_coll.cpp.
|
inlinevirtual |
Reimplemented in mfem::L2_FECollection.
Definition at line 67 of file fe_coll.hpp.
1.8.5