In this page, we detail the different hexahedral finite elements currently implemented in Montjoie. The reference element is the unit cube, look at the figure below for the conventions used

Unit cube

If you want to implement a new hexahedral element, you have to derive the class HexahedronReference. This class is the base class for all hexahedral elements implemented in Montjoie. It depends on the type of element (1: scalar finite element, 2: edge elements, 3: facet elements).

// for a new hexahedral finite element (1 : scalar)
class MyFiniteElement : public HexahedronReference<1>
{
   public :
   
   // number of dofs on each edge, quadrangle, hexahedron, etc ?
   void ConstructNumberMap(NumberMap& map, int dg) const;

   // main method used to construct the finite element object
   // r : order of approximation for basis functions
   // rgeom : order of approximation for geometry (shape functions)
   // rquad : order of quadrature
   // if rgeom is equal to 0 (and similarly rquad equal to 0)
   // we are considering that rgeom = r, and rquad = r
   // if type_quad is equal to -1, we use the default quadrature
   // formulas (Gauss-Legendre formulas)
   // rsurf_tri, rsurf_quad : orders for boundary integrals
   // type_surf_tri, type_surf_quad : type of integration for triangles and quadrangles
   void ConstructFiniteElement(int r, int rgeom = 0, int rquad = 0, int type_quad = -1,
                               int rsurf_tri = 0, int rsurf_quad = 0, int type_surf_tri = -1,
                               int type_surf_quad = -1, int gauss_z = -1);

   // computation of phi_i(x) for a given point x
   // the values phi_i(x) have to be placed in the array phi
   void ComputeValuesPhiRef(const R3& x, VectReal_wp& phi) const;

   // computation of  \nabla phi_i(x) for a given point x
   void ComputeGradientPhiRef(const R3& x, VectR3& grad_phi) const;

};

HexahedronGeomReference

The class HexahedronGeomReference implements shape functions methods for hexahedron elements. These shape functions are computed by tensorization of 1-D Lagrange polynomials :

where xk are 1-D nodal points. These points are equal to Gauss-Lobatto points. Most of the methods of this class come from the parent class ElementGeomReference. Below, we detail the methods and attributes specific to the class HexahedronGeomReference.

GetNbPointsNodalInside returns the number of nodal points inside the hexahedron.

HexahedronReference

The class HexahedronReference is the base class for all hexahedral finite elements. The methods come from the parent class ElementReference.

HexahedronGauss

This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The considered finite element space consists of polynomials of degree less or equal to r for each variable :

A tensorized basis is used

with Lagrange interpolation function :

The interpolation points are chosen as Gauss-Lobatto points. Gauss-Legendre points are used for quadrature formulas (both for volume integrals and boundary integrals). The class HexahedronGauss allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidClassical, WedgeClassical, TetrahedronClassical. The class HexahedronGauss is the class used for hexahedra if you have TypeElement = TETRAHEDRON_CLASSICAL in the datafile.

HexahedronLobatto

This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The basis functions are based on Gauss-Lobatto points, they are the same as for HexahedronGauss. The quadrature formulas are based on Gauss-Lobatto points (instead of Gauss-Legendre for HexahedronGauss), leading to a diagonal mass matrix (mass lumping). Boundary integrals are also computed with Gauss-Lobatto points. Finite elements computations are faster than for HexahedronGauss, and elementary matrices sparser. This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidLobatto, WedgeLobatto, TetrahedronClassical. The class HexahedronLobatto is the class used for hexahedra if you have TypeElement = TETRAHEDRON_LOBATTO or TypeElement = HEXAHEDRON_LOBATTO in the datafile.

HexahedronDgGauss

This class implements a scalar finite element with orthogonal basis functions. It can be used to discretize L2 (with discontinuous Galerkin formulations). The basis functions are based on Gauss-Legendre points instead of Gauss-Lobatto points, which is possible since this element is used with a discontinuous Galerkin method. Quadrature formulas are also based on Gauss points, leading to a diagonal mass matrix. This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidClassical (or PyramidDgOrtho), WedgeClassical (or WedgeDgOrtho), TetrahedronClassical (or TetrahedronDgOrtho). The class HexahedronDgGauss is the class used for hexahedra if you have TypeElement = TETRAHEDRON_DG_CLASSICAL or TypeElement = TETRAHEDRON_DG_ORTHO in the datafile.

HexahedronDgLegendre

This class implements a scalar finite element with orthogonal basis functions. It can be used to discretize L2 (with discontinuous Galerkin formulations). The basis functions used for this class are the orthogonal functions (generating Pr instead of Qr) :

Because of tensorized structure, fast computations are obtained for a large order of approximation. No quadrature rules are used for evaluating integrals, a quadrature-free approach is used to obtain faster algorithms. This class can be used only if the mesh contains cubes (and not general hexahedra) in order to use this quadrature-free approach. This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidDgLegendre, WedgeDgLegendre, TetrahedronDgOrtho. The class HexahedronDgLegendre is the class used for hexahedra if you have TypeElement = TETRAHEDRON_DG_LEGENDRE in the datafile.

HexahedronHierarchic

This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The finite element is Qr (same as for HexahedronGauss). The basis functions for this class are hierarchical, allowing the use of a variable order in a mesh. The basis functions are equal to :

The tensor structure of this basis is exploited to obtain faster computations for large values of r. Gauss quadrature rules are used for volume integrals and boundary integrals. The class HexahedronHierarchic allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHierarchic, WedgeHierarchic, TetrahedronHierarchic. The class HexahedronHierarchic is the class used for hexahedra if you have TypeElement = TETRAHEDRON_HIERARCHIC in the datafile.

HexahedronHcurlFirstFamily

This class implements interpolatory basis functions for edge elements of Nedelec's first family. It can be used to discretize H(curl) space. The finite element space considered here is equal to

The basis functions are based on r Gauss points and r+1 Gauss-Lobatto points :

The basis functions are then equal to :

Gauss-Lobatto points are used for the quadrature formulas, leading to fast computation of the elementary matrices. The class HexahedronHcurlFirstFamily allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHcurlFirstFamily, WedgeHcurlFirstFamily, TetrahedronHcurlFirstFamily. The class HexahedronHcurlFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_FIRST_FAMILY or TypeElement = HEXAHEDRON_FIRST_FAMILY in the datafile.

HexahedronHcurlOptimalFirstFamily

This class implements interpolatory basis functions for optimal hexahedral edge elements of Nedelec's first family. It can be used to discretize H(curl) space. The involved finite element space is equal to

With this finite element space, the convergence towards the solution is optimal in O(hr) in H(curl) norm. If you use HexahedronHcurlFirstFamily, the convergence for general meshes will be in O(hr-1). The basis functions are based on r+2 Gauss-Lobatto points and the r interior Gauss-Lobatto points :

The basis functions are then equal to :

The quadrature formulas are based on r+2 Gauss-Lobatto points. This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHcurlOptimalFirstFamily, WedgeHcurlOptimalFirstFamily, TetrahedronHcurlOptimalFirstFamily. The class HexahedronHcurlOptimalFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_OPTIMAL_FIRST_FAMILY in the datafile.

HexahedronHcurlLobatto

This finite element class implements interpolatory basis functions for hexahedral edge elements of Nedelec's second family. It can be used to discretize H(curl) space. The finite element space is given as

Gauss-Lobatto points are used both for basis functions and quadrature rules (leading to a block-diagonal mass matrix and efficient computations). However this type of element produces spurious modes, particularly for distorted meshes. Numerical computations cannot be conducted over hybrid meshes, however you can use purely tetrahedral or purely hexahedral meshes, the corresponding class for tetrahedra is the class TetrahedronHcurlLobatto.

HexahedronHcurlHpFirstFamily

This class implements hierarchical basis functions for edge hexahedral edge elements of Nedelec's first family. It can be used to discretize H(curl) space. The finite element space is the same as for HexahedronHcurlFirstFamily. By using hierarchical functions, a variable order can be set for each element of the mesh. The basis functions are equal to :

This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHcurlHpFirstFamily, WedgeHcurlHpFirstFamily, TetrahedronHcurlOptimalHpFirstFamily. The class HexahedronHcurlHpFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_HP_FIRST_FAMILY in the datafile.

HexahedronHcurlOptimalHpFirstFamily

This class implements hierarchical basis functions for edge hexahedral optimal edge elements of Nedelec's first family. It can be used to discretize H(curl) space. The finite element space is the same as for HexahedronHcurlOptimalFirstFamily. The basis functions are equal to :

This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHcurlOptimalHpFirstFamily, WedgeHcurlOptimalHpFirstFamily, TetrahedronHcurlOptimalHpFirstFamily. The class HexahedronHcurlOpitmalHpFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_OPTIMAL_HP_FIRST_FAMILY in the datafile.

HexahedronHdivFirstFamily

This class implements interpolatory basis functions for facet elements of Nedelec's first family. It can be used to discretize H(div) space. The finite element space considered here is equal to

The basis functions are based on r Gauss points and r+1 Gauss-Lobatto points :

The basis functions are then equal to :

Gauss-Lobatto points are used for the quadrature formulas, leading to fast computation of the elementary matrices. The class HexahedronHdivFirstFamily allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHdivFirstFamily, WedgeHdivFirstFamily, TetrahedronHdivFirstFamily. The class HexahedronHdivFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_FIRST_FAMILY in the datafile.

HexahedronHdivOptimalFirstFamily

This class implements interpolatory basis functions for facet elements of optimal Nedelec's first family. It can be used to discretize H(div) space. The finite element space considered here is equal to

With this finite element space, the convergence towards the solution is optimal in O(hr) in H(div) norm. If you use HexahedronHdivFirstFamily, the convergence for general meshes will be in O(hr-2). The basis functions are based on r+1 Gauss points and r+3 Gauss-Lobatto points :

The basis functions are then equal to :

Gauss-Lobatto points are used for the quadrature formulas, leading to fast computation of the elementary matrices. The class HexahedronHdivFirstFamily allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHdivOptimalFirstFamily, WedgeHdivOptimalFirstFamily, TetrahedronHdivFirstFamily. The class HexahedronHdivOptimalFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_OPTIMAL_FIRST_FAMILY in the datafile.

HexahedronHdivHpFirstFamily

This class implements hierarchical basis functions for hexahedral facet elements of Nedelec's first family. It can be used to discretize H(div) space. The finite element space is the same as for HexahedronHdivFirstFamily. By using hierarchical functions, a variable order can be set for each element of the mesh. The basis functions are equal to :

This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHdivHpFirstFamily, WedgeHdivHpFirstFamily, TetrahedronHdivOptimalHpFirstFamily. The class HexahedronHdivHpFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_HP_FIRST_FAMILY in the datafile.

HexahedronHdivOptimalHpFirstFamily

This class implements hierarchical basis functions for hexahedral facet elements of optimal Nedelec's first family. It can be used to discretize H(div) space. The finite element space is the same as for HexahedronHdivOptimalFirstFamily. By using hierarchical functions, a variable order can be set for each element of the mesh. The basis functions are equal to :

This class allows numerical computations over hybrid meshes, the corresponding classes for pyramids, prisms and tetrahedra are PyramidHdivOptimalHpFirstFamily, WedgeHdivOptimalHpFirstFamily, TetrahedronHdivOptimalHpFirstFamily. The class HexahedronHdivOptimalHpFirstFamily is the class used for hexahedra if you have TypeElement = TETRAHEDRON_OPTIMAL_HP_FIRST_FAMILY in the datafile.

GetNbPointsNodalInside

Syntax

int GetNbPointsNodalInside() const

This method returns the number of nodal points strictly inside the unit cube. It should be equal to (r-1)3 where r is the polynomial degree.

Example :

int r = 4; // polynomial degree
HexahedronGeomReference hex;
hex.ConstructFiniteElement(r);

// number of nodal points inside the hexahedron ?
int nb_nodes_inside = hex.GetNbPointsNodalInside(); // should give 27 here

Location :

FiniteElement/Hexahedron/HexahedronGeomReference.hxx   FiniteElement/Hexahedron/HexahedronGeomReferenceInline.cxx