Quadrangular finite elements spaces
In this page, we detail the different quadrangular finite elements currently implemented in Montjoie. The reference element is the unit square, look at the figure below for the conventions used
If you want to implement a new quadrangular element, you have to derive the class QuadrangleReference. This class is the base class for all triangular elements implemented in Montjoie. It depends on the type of element (1: scalar finite element, 2: edge elements, 3: facet elements).
// for a new triangular finite element (1 : scalar) class MyFiniteElement : public QuadrangleReference<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 formulas for quadrilaterals) // redge, type_edge : quadrature formulas for edges void ConstructFiniteElement(int r, int rgeom = 0, int rquad = 0, int type_quad = -1, int redge = 0, int type_edge = -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 R2& x, VectReal_wp& phi) const; // computation of \nabla phi_i(x) for a given point x void ComputeGradientPhiRef(const R2& x, VectR2& grad_phi) const; };
QuadrangleGeomReference
The class QuadrangleGeomReference implements shape functions methods for quadrangular 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 QuadrangleGeomReference.
GetNbPointsNodalInside | returns the number of nodal points inside the quadrangle. |
SetNodalPoints | constructs shapes functions from a set of nodal points |
Methods and attributes specific to QuadrangleReference (inherited from ElementReference)
The class QuadrangleReference is the base class for all quadrangular finite elements. Below we list methods and attributes specific to this class (other methods come from the parent class ElementReference).
ConstructQuadrature | constructs the quadrature rules for 2-D and boundary integrals. |
ConstructFiniteElement1D | constructs the boundary finite element (trace on an edge) |
QuadrangleGauss
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. In the picture below, we have displayed the degrees of freedom for r=3, 6 and 9:
The used quadrature formulas are Gauss-Legendre for 2-D integrals and for boundary integrals. The class QuadrangleGauss allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleClassical. The class QuadrangleGauss is the class used for quadrilaterals if you have TypeElement = TRIANGLE_CLASSICAL in the datafile.
QuadrangleLobatto
This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The basis functions are the same for QuadrangleGauss and QuadrangleLobatto. The only difference between the two classes is that integrals are evaluated with Gauss-Lobatto points with QuadrangleLobatto. This leads to mass lumping (the mass matrix is diagonal) and faster computations. Besides, the class QuadrangleLobatto allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleLobatto. The class QuadrangleLobatto is the class used for quadrilaterals if you have TypeElement = TRIANGLE_LOBATTO or TypeElement = QUADRANGLE_LOBATTO in the datafile.
QuadrangleDgGauss
This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The finite element space is Qr (identic to QuadrangleGauss). This class implements nodal basis functions on quadrilateral elements for a Discontinuous Galerkin formulation. The basis functions are similar to basis functions of QuadrangleGauss, but they are based on Gauss points instead of Gauss-Lobatto points. The integrals are also evaluated with Gauss points, leading to a diagonal mass matrix. Besides, the class QuadrangleDgGauss allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleClassical. The class QuadrangleDgGauss is the class used for quadrilaterals if you have TypeElement = QUADRANGLE_DG_GAUSS in the datafile.
QuadrangleHierarchic
This class implements a scalar finite element with interpolatory basis functions. It can be used to discretize H1. The finite element space is Qr (identic to QuadrangleGauss). This class implements hierarchical basis functions on quadrilateral elements. The basis functions are the following :
where Pi1, 1 are Jacobi polynomials orthogonal with respect to the weight (1-x)(1+x). The class QuadrangleHierarchic allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHierarchic. The class QuadrangleHierarchic is the class used for quadrangles if you have TypeElement = TRIANGLE_HIERARCHIC in the datafile.
QuadrangleHcurlFirstFamily
This class implements interpolatory basis functions for quadrilateral edge elements of Nedelec's first family. It can be used to discretize H(curl). The involved finite element space 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 :
In the picture below, degrees of freedom are represented with their orientation for r=2, 4 and 6.
The class QuadrangleHcurlFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHcurlFirstFamily. The class QuadrangleHcurlFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_FIRST_FAMILY or TypeElement = QUADRANGLE_FIRST_FAMILY in the datafile.
QuadrangleHcurlOptimalFirstFamily
This class implements interpolatory basis functions for optimal quadrilateral edge elements of Nedelec's first family. It can be used to discretize H(curl). 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 QuadrangleHcurlFirstFamily, 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 integrals are computed with the (r+2)2 Gauss-Lobatto points. The class QuadrangleHcurlOptimalFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHcurlFirstFamily. The class QuadrangleHcurlOptimalFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_OPTIMAL_FIRST_FAMILY in the datafile.
QuadrangleHcurlLobatto
This class implements interpolatory basis functions for quadrilateral edge elements of Nedelec's second family. It can be used to discretize H(curl). The involved finite element space is equal to
The basis functions are similar to the ones of QuadrangleLobatto with an orientation (along x or y). It should be noticed that the use of Nedelec's second family on quadrilateral elements provides spurious modes, especially on highly distorted meshes. The class QuadrangleHcurlLobatto allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHcurlSecondFamily. The class QuadrangleHcurlLobatto is the class used for quadrangles if you have TypeElement = QUADRANGLE_HCURL_LOBATTO in the datafile.
QuadrangleHcurlHpFirstFamily
This class implements hierarchical basis functions for quadrilateral edge elements of Nedelec's first family. It can be used to discretize H(curl). The involved finite element space is equal to
The basis functions are equal to :
The class QuadrangleHcurlHpFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHcurlOptimalHpFirstFamily. The class QuadrangleHcurlHpFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_HP_FIRST_FAMILY in the data file.
QuadrangleHcurlOptimalHpFirstFamily
This class implements hierarchical basis functions for optimal quadrilateral edge elements of Nedelec's first family. It can be used to discretize H(curl). The involved finite element space is equal to
The basis functions are equal to
The class QuadrangleHcurlOptimalHpFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHcurlOptimalHpFirstFamily. The class QuadrangleHcurlOptimalHpFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_OPTIMAL_HP_FIRST_FAMILY in the data file.
QuadrangleHdivFirstFamily
This class implements interpolatory basis functions for quadrilateral facet elements of Nedelec's first family. It can be used to discretize H(div). The involved finite element space 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 :
In the picture below, degrees of freedom are represented with their orientation for r=2, 4 and 6.
The class QuadrangleHdivFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHdivFirstFamily. The class QuadrangleHdivFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_FIRST_FAMILY or TypeElement = QUADRANGLE_FIRST_FAMILY in the datafile.
QuadrangleHdivOptimalFirstFamily
This class implements interpolatory basis functions for optimal quadrilateral facet elements of Nedelec's first family. It can be used to discretize H(div). 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(div) norm. If you use QuadrangleHdivFirstFamily, 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 integrals are computed with the (r+2)2 Gauss-Lobatto points. The class QuadrangleHdivOptimalFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHdivFirstFamily. The class QuadrangleHdivOptimalFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_OPTIMAL_FIRST_FAMILY in the datafile.
QuadrangleHdivHpFirstFamily
This class implements hierarchical basis functions for quadrilateral facet elements of Nedelec's first family. It can be used to discretize H(div). The involved finite element space is equal to
The basis functions are equal to :
The class QuadrangleHdivHpFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHdivOptimalHpFirstFamily. The class QuadrangleHdivHpFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_HP_FIRST_FAMILY in the data file.
QuadrangleHdivOptimalHpFirstFamily
This class implements hierarchical basis functions for optimal quadrilateral facet elements of Nedelec's first family. It can be used to discretize H(div). The involved finite element space is equal to
The basis functions are equal to
The class QuadrangleHdivOptimalHpFirstFamily allows numerical computations over hybrid meshes, the corresponding class for triangles is TriangleHdivOptimalHpFirstFamily. The class QuadrangleHdivOptimalHpFirstFamily is the class used for quadrangles if you have TypeElement = TRIANGLE_OPTIMAL_HP_FIRST_FAMILY in the data file.
GetNbPointsNodalInside
Syntax
int GetNbPointsNodalInside() const |
This method returns the number of nodal points strictly inside the unit square. It should be equal to (r-1)(r-1) where r is the polynomial degree.
Example :
int r = 4; // polynomial degree QuadrangleGeomReference quad; quad.ConstructFiniteElement(r); // number of nodal points inside the quadrangle ? int nb_nodes_inside = quad.GetNbPointsNodalInside(); // should give 9 here
Location :
FiniteElement/Quadrangle/QuadrangleGeomReference.hxx FiniteElement/Quadrangle/QuadrangleGeomReferenceInline.cxx
SetNodalPoints
Syntax
void SetNodalPoints(int r, const VectR2& pts) const |
This method constructs the shape function from nodal points provided by the user. It can be used in replacement of ConstructFiniteElement if you want to use your own set of points (instead of Gauss-Lobatto points). The number of points should be equal to (r+1)(r+1) with r+1 nodal points on each edge.
Example :
int r = 3; VectR2 pts((r+1)*(r+1)); pts(0).Init(0.0, 0.0); pts(1).Init(1.0, 0.0); pts(2).Init(1.0, 1.0); pts(3).Init(0.0, 1.0); // etc QuadrangleGeomReference quad; // you can provide your own set of nodal points quad.SetNodalPoints(r, pts); // evaluating shape functions R2 pt_loc(0.2, 0.4); VectReal_wp phi; quad.ComputeValuesPhiNodalRef(pt_loc, phi);
Location :
FiniteElement/Quadrangle/QuadrangleGeomReference.hxx FiniteElement/Quadrangle/QuadrangleGeomReference.cxx
ConstructQuadrature
Syntax
void ConstructQuadrature( | int r, int type_quad, int r_edge) |
This method computes the quadrature rules to use for the interior and for the boundary. This method should not be called in a regular use, since it is already called by ConstructFiniteElement. If you define a new finite element that derives from QuadrangleReference, you can call this method to construct quadrature rules for the quadrangle. If r is the order for quadrature rules, the rules should integrate exactly polynomials of degree 2r+1 or 2r-1 such that the elementary matrices are integrated accurately.
Parameters
- r (in)
- order for quadrature rules for the quadrangle
- type_quad (in)
- type of integration rule as defined by the method Globatto::ConstructQuadrature
- r_edge (in)
- order for quadrature rules for the boundary
Example :
class MyQuadrangle : public QuadrangleReference<1> { public: void ConstructFiniteElement(int r, int rgeom = 0, int rquad = 0, int type_quad = -1, int rsurf = 0, int type_surf = -1) { // constructing quadrature rules : this->ConstructQuadrature(r, type_quad, rsurf); } };
Location :
FiniteElement/Quadrangle/QuadrangleReference.hxx FiniteElement/Quadrangle/QuadrangleReference.cxx
ConstructFiniteElement1D
Syntax
void ConstructFiniteElement1D( | int r, int rgeom, int rquad, int type_quad) |
This virtual method computes the 1-D finite element which is the trace of the quadrangular finite element. The default constructed 1-D element is EdgeGaussReference, which corresponds to the trace of QuadrangleGauss. Since this method is virtual, you can overload it without overloading ConstructFiniteElement.
Parameters
- r (in)
- order of approximation
- rgeom (in)
- order for geometry
- rquad (in)
- order for quadrature rules for the unit interval
- type_quad (in)
- type of rule for the edges as defined by the method Globatto::ConstructQuadrature
Example :
class MyQuadrangle : public QuadrangleReference<1> { public: void ConstructFiniteElement1D(int r, int rgeom, int rquad, int type_quad) { // constructing 1-D finite element, e.g. EdgeLobattoReference if (this->element_surface != NULL) delete this->element_surface; EdgeLobattoReference* edge = new EdgeLobattoReference(); edge->ConstructFiniteElement(r, rgeom, rquad, type_quad); this->element_surface = edge; } };
Location :
FiniteElement/Quadrangle/QuadrangleReference.hxx FiniteElement/Quadrangle/QuadrangleReference.cxx