12 #ifndef MFEM_FE_COLLECTION
13 #define MFEM_FE_COLLECTION
15 #include "../config/config.hpp"
29 template <Geometry::Type geom>
30 static inline void GetNVE(
int &nv,
int &ne);
32 template <Geometry::Type geom,
typename v_t>
33 static inline void GetEdge(
int &nv, v_t &v,
int &ne,
int &e,
int &eo,
37 typename v_t,
typename e_t,
typename eo_t>
38 static inline void GetFace(
int &nv, v_t &v,
int &ne, e_t &e, eo_t &eo,
39 int &nf,
int &f,
int &fg,
int &fo,
50 virtual const char *
Name()
const {
return "Undefined"; }
96 {
return H1_dof[GeomType]; }
103 const int *
GetDofMap(
int GeomType)
const;
145 {
return L2_Elements[GeomType]; }
148 if (L2_Elements[GeomType])
150 return L2_Elements[GeomType]->
GetDof();
155 virtual const char *
Name()
const {
return d_name; }
160 return Tr_Elements[GeomType];
182 void InitFaces(
const int p,
const int dim,
const int map_type,
199 {
return RT_dof[GeomType]; }
246 {
return ND_dof[GeomType]; }
279 mutable char name[16];
291 QuadrilateralFE->
Reset();
292 ParallelepipedFE->
Reset();
310 virtual const char *
Name()
const {
return name; }
338 virtual const char *
Name()
const {
return "Linear"; }
362 virtual const char *
Name()
const {
return "Quadratic"; }
382 virtual const char *
Name()
const {
return "QuadraticPos"; }
406 virtual const char *
Name()
const {
return "Cubic"; }
426 virtual const char *
Name()
const {
return "CrouzeixRaviart"; }
448 virtual const char *
Name()
const {
return "LinearNonConf3D"; }
470 virtual const char *
Name()
const {
return "RT0_2D"; }
491 virtual const char *
Name()
const {
return "RT1_2D"; }
512 virtual const char *
Name()
const {
return "RT2_2D"; }
532 virtual const char *
Name()
const {
return "Const2D"; }
553 virtual const char *
Name()
const {
return "LinearDiscont2D"; }
574 virtual const char *
Name()
const {
return "GaussLinearDiscont2D"; }
588 virtual const char *
Name()
const {
return "P1OnQuad"; }
609 virtual const char *
Name()
const {
return "QuadraticDiscont2D"; }
625 virtual const char *
Name()
const {
return "QuadraticPosDiscont2D"; }
646 virtual const char *
Name()
const {
return "GaussQuadraticDiscont2D"; }
667 virtual const char *
Name()
const {
return "CubicDiscont2D"; }
688 virtual const char *
Name()
const {
return "Const3D"; }
709 virtual const char *
Name()
const {
return "LinearDiscont3D"; }
730 virtual const char *
Name()
const {
return "QuadraticDiscont3D"; }
754 virtual const char *
Name()
const {
return "RefinedLinear"; }
775 virtual const char *
Name()
const {
return "ND1_3D"; }
797 virtual const char *
Name()
const {
return "RT0_3D"; }
818 virtual const char *
Name()
const {
return "RT1_3D"; }
833 {
return (GeomType == _GeomType) ? Local_Element : NULL; }
835 {
return (GeomType == _GeomType) ? Local_Element->
GetDof() : 0; }
838 virtual const char *
Name()
const {
return d_name; }
Abstract class for Finite Elements.
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
virtual const char * Name() const
virtual const char * Name() const
virtual const char * Name() const
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
For scalar fields; preserves point values.
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
int RT_dof[Geometry::NumGeom]
FiniteElement * ND_Elements[Geometry::NumGeom]
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
virtual const char * Name() const
virtual FiniteElementCollection * GetTraceCollection() const
virtual const char * Name() const
GaussQuadraticDiscont2DFECollection()
virtual ~NURBSFECollection()
void SubDofOrder(int Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
RT_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
FiniteElementCollection * GetTraceCollection() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const char * Name() const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
QuadraticDiscont3DFECollection()
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual ~L2_FECollection()
LinearNonConf3DFECollection()
virtual const char * Name() const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
RefinedLinearFECollection()
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Piecewise-(bi)linear continuous finite elements.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Class for quadratic FE on triangle.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Class for quadratic FE on interval.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
FiniteElementCollection * GetTraceCollection() const
virtual int DofForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
FiniteElementCollection * GetTraceCollection() const
virtual const char * Name() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE, const int ob_type=BasisType::GaussLegendre)
virtual const char * Name() const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual const char * Name() const
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Finite element collection on a macro-element.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Class for refined bi-linear FE on quadrilateral.
Possible basis types. Note that not all elements can use all BasisType(s).
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
QuadraticPosDiscont2DFECollection()
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
virtual const char * Name() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual ~RT_FECollection()
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
virtual int DofForGeometry(int _GeomType) const
virtual int DofForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Piecewise-(bi)cubic continuous finite elements.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const char * Name() const
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Class for refined linear FE on interval.
virtual int DofForGeometry(int GeomType) const
int H1_dof[Geometry::NumGeom]
Class for refined linear FE on triangle.
Class for refined linear FE on tetrahedron.
FiniteElement * RT_Elements[Geometry::NumGeom]
Class for constant FE on triangle.
static void GetNVE(int &nv, int &ne)
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
For scalar fields; preserves volume integrals.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
ND_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Class for refined trilinear FE on a hexahedron.
Discontinuous collection defined locally by a given finite element.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
CubicDiscont2DFECollection()
Version of QuadraticDiscont2DFECollection with positive basis functions.
virtual ~ND_FECollection()
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Class for quadratic FE on tetrahedron.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
QuadraticDiscont2DFECollection()
const int * GetDofMap(int GeomType) const
Get the Cartesian to local H1 dof map.
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Version of QuadraticFECollection with positive basis functions.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual ~Local_FECollection()
virtual const char * Name() const
virtual int DofForGeometry(int GeomType) const
virtual ~H1_FECollection()
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
LinearDiscont3DFECollection()
Class for linear FE on tetrahedron.
Class for linear FE on interval.
Crouzeix-Raviart finite element on quadrilateral.
Class for linear FE on triangle.
LinearDiscont2DFECollection()
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Class for quadratic FE on triangle with nodes at the "Gaussian" points.
virtual const char * Name() const
virtual const char * Name() const
FiniteElementCollection * GetTraceCollection() const
FiniteElement * H1_Elements[Geometry::NumGeom]
Piecewise-linear nonconforming finite elements in 3D.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Crouzeix-Raviart nonconforming elements in 2D.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
H1Pos_FECollection(const int p, const int dim=3)
virtual int DofForGeometry(int GeomType) const
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
virtual const char * Name() const
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
void SetOrder(int Order) const
Set the order and the name, based on the given Order: either a positive number for fixed order...
virtual const char * Name() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const char * Name() const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Class for bi-quadratic FE on quadrilateral.
Class for tri-linear FE on cube.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual const char * Name() const
Class for bilinear FE on quadrilateral.
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual const char * Name() const
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Class for linear FE on triangle with nodes at the 3 "Gaussian" points.
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
GaussLinearDiscont2DFECollection()
virtual int DofForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
virtual const char * Name() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
QuadraticPosFECollection()
Tensor products of 1D FEs (only degree 2 is functional)
Piecewise-(bi)quadratic continuous finite elements.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
int HasFaceDofs(int GeomType) const
Arbitrary order H(curl)-conforming Nedelec finite elements.
Class for cubic FE on tetrahedron.
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, int &fg, int &fo, const int face_info)
virtual const char * Name() const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Arbitrary order H1-conforming (continuous) finite elements.
virtual const char * Name() const
virtual const char * Name() const
virtual const FiniteElement * FiniteElementForGeometry(int _GeomType) const
CrouzeixRaviartFECollection()
virtual int DofForGeometry(int GeomType) const =0
virtual int * DofOrderForOrientation(int GeomType, int Or) const
virtual int DofForGeometry(int GeomType) const
Local_FECollection(const char *fe_name)
virtual int DofForGeometry(int GeomType) const
virtual const char * Name() const
virtual int DofForGeometry(int GeomType) const
Linear (P1) finite elements on quadrilaterals.
virtual int DofForGeometry(int GeomType) const
virtual ~FiniteElementCollection()
int ND_dof[Geometry::NumGeom]
virtual int DofForGeometry(int GeomType) const
Crouzeix-Raviart finite element on triangle.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
Arbitrary order "L2-conforming" discontinuous finite elements.
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const