Appendix A. Finite element method list
Let us recall that all finite element methods defined in
getfem++ are declared in the file getfem_fem.h and that a descriptor on a finite
element method is obtained thanks to the function
getfem::pfem pf =
getfem::fem_descriptor("name of method");
where "name of method" is a
string to be choosen among the existing methods.
Symbols representing degree of freedom types
It is possible to define a classical PK
Lagrange element of arbitrary dimension and arbitrary degree. Each
degree of freedom of such an element corresponds to the value of
the function on a corresponding node. The grid of node is the
so-called Lagrange grid. Figures *,
* and * show examples in dimension 1, 2 and
3.

Examples of classical PK Lagrange
elements on a segment.
 |
 |
| P1 element, 3
d.o.f., C0 |
P2 element, 6
d.o.f., C0 |
|
 |
 |
| P3 element, 10
d.o.f., C0 |
P6 element, 28
d.o.f., C0 |
Examples of classical PK Lagrange
elements on a triangle.
The number of degree of freedom for a classical
PK Lagrange element of dimension P and
degree K is ((P+K)!)/(P!K!). For instance, in
dimension 2 (P = 2), this value is ((P+1) (P+2))/(2)
and in dimension 3 (P = 3), it is ((P+1) (P+2)
(P+3))/(6).
 |
 |
| P1 element, 4
d.o.f., C0, |
P2 element, 10
d.o.f., C0 |
 |
|
| P4 element, 35
d.o.f., C0 |
|
Examples of classical PK Lagrange
elements on a tetrahedron.
The particular way used in getfem++ to numerate the nodes
are also shown in figures *, * and *.
Using another numeration, let
be somme indices such that
0 <= i0, i1, ...
iP <= K, and ∑n =
0P in = K.
Then, the coordinate of a node can be computed as
ai0, i1, ...
iP = ∑n = 0P
(in)/(K)Sn, for K != 0,
where
S0, S1, ... SN are
the vertices of the simplex (for
K = 0 the particular choice
a0, 0, ... 0 = ∑n = 0P
(1)/(P+1)Sn has been chosen). Then each base
function, corresponding of each node
ai0,
i1, ... iP is defined by
φi0, i1, ...
iP = Πn = 0P
Πj=0in-1 ((K
λn - j)/(j+1)).
where
λn are the barycentric coordinates,
i.e. the polynomials of degree 1 whose value is
1 on the
vertex
Sn and whose value is
0 on other
vertices. On the reference element, one has
λP = 1 - x0 -
x1 - ... - xP-1.
When between two elements of the same degrees (even with different
dimensions), the d.o.f. of a common face are linked, the element is
of class
C0. This means that the global
polynomial is continuous. If you try to link elements of different
degrees, you will get some trouble with the unlinked d.o.f. This is
not automatically supported by getfem++ , so you will have to
support it (add constraints on these d.o.f.).
For some applications (computation of a gradient for instance) one
may not want the d.o.f. of a common face to be linked. This is why
there are two versions of the classical
PK
Lagrange element.
| Classical PK
Lagrange element |
| "FEM_PK(P, K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K, 0
<= K <= 255 |
P,
1 <= P <=
255 |
((K+P)!)/(K! P!) |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
| Discontinuous
PK Lagrange element |
| "FEM_PK_DISCONTINUOUS(P,
K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K, 0
<= K <= 255 |
P,
1 <= P <=
255 |
((K+P)!)/(K! P!) |
discontinuous |
No (Q = 1) |
Yes (M = Id) |
Yes |
Even thought Lagrange elements are defined for arbitrary degrees,
to choose a hight degree can be problematic for a large number of
applications due to the "noisy" caracteristic of the lagrange
basis. Those element are recommended for the basic interpolation
but for p.d.e. applications elements with hierarchical basis are
preferable (see the corresponding section).
Classical Lagrange elements on parallelepipeds or prisms are
obtained as tensor product of Lagrange elements on simplices. When
two element are defined, one on a dimension P1
and the other in dimension P2, one obtains the
base functions of the tensorial product (on the reference element)
as
φ'ij(x,y) =
φ'1i(x)
φ'2j(y), x in
ℜP1, y in
ℜP2,
where
φ'1i and
φ'2i are respectively the base
functions of the first and second element.
 |
 |
| Q1 element, 4
d.o.f., C0 |
Q3 element, 16
d.o.f., C0 |
|
Examples of classical QK Lagrange
elements in dimension 2
The QK element on a parallelepiped of
dimension P is obtained as the tensorial product of P
classical PK element on the segment. Examples in
dimension 2 are shown in figure * and in dimension 3 in figure
*.
A prism in dimension P > 1 is the direct product of a
simplex of dimension P-1 with a segment. The
PK ⊗PK element on this prism is
the tensorial product of the classical PK element
on a simplex of dimension P-1 with the classical
PK element on a segment. For P=2 this
coincide with a parallelepiped. Examples in dimension 3 are
shown in figure *. This is also
possible not to have the same degree on each dimension. An example
is shown on figure *.
 |
 |
| Q1 element, 8
d.o.f., C0 |
Q3 element, 64
d.o.f., C0 |
 |
 |
| P1
⊗P1 element, 6 d.o.f.,
C0 |
P3
⊗P3 element, 40 d.o.f.,
C0 |
|
Examples of classical Lagrange elements in
dimension 3
P2 ⊗P1 Lagrange
element on a prism, 12 d.o.f., C0
| QK Lagrange
element on parallelepipeds |
| "FEM_QK(P, K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| KP, 0
<= K <= 255 |
P,
2 <= P <=
255 |
(K+1)P |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
| PK
⊗PK Lagrange element on prisms |
| "FEM_PK_PRISM(P, K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 2K, 0
<= K <= 255 |
P,
2 <= P <=
255 |
(K+1)
× ((K+P-1)!)/(K! (P-1)!) |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
| PK1
⊗PK2 Lagrange element on
prisms |
| "FEM_PRODUCT(FEM_PK(P-1,
K1), FEM_PK(1, K2))" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K1+K2, 0 <= K1,K2 <=
255 |
P,
2 <= P <=
255 |
(K2+1)
× ((K1+P-1)!)/(K1! (P-1)!) |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |

Incomplete Q2 element, 8 d.o.f.,
C0
| Incomplete Q2
Lagrange element on quadrilateral (Quad 8 serendipity
element) |
| "FEM_INCOMPLETE_Q2" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
8 |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
The idea behind hierarchical basis is the description of the
solution at different level : a rought level, a more refined level
... In the same discretisation some degrees of freedom represent
the rought description, some other the more rafined and so on. This
correspond to imbricated spaces of discretisation. The hierarchical
basis contains a basis of each of these spaces (this is not the
case in classical Lagrange elements when the mesh is
refined).
Among the advantages, the condition number of rigidity matrices can
be greatly improved, it allows local raffinement and a resolution
with a multigrid approach.
PK Hierarchical element on a segment,
C0
| PK Classical
Lagrange element on simplices but with a hierarchical basis with
respect to the degree |
|
"FEM_PK_HIERARCHICAL(P,K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K,
0 <= K <=
255 |
P,
1 <= P <=
255 |
((K+P)!)/(K! P!) |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
| QK Classical
Lagrange element on parallelepipeds but with a hierarchical basis
with respect to the degree |
|
"FEM_QK_HIERARCHICAL(P,K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K,
0 <= K <=
255 |
P,
2 <= P <=
255 |
(K+1)P |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
| PK Classical
Lagrange element on prisms but with a hierarchical basis with
respect to the degree |
|
"FEM_PK_PRISM_HIERARCHICAL(P,K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K,
0 <= K <=
255 |
P,
2 <= P <=
255 |
(K+1)
× ((K+P-1)!)/(K! (P-1)!) |
C0 |
No (Q = 1) |
Yes (M = Id) |
Yes |
some particular choices :
P4 will be build with
the basis of the
P1, the additional basis of the
P2 then the additionnal basis of the
P4.
P6 will be build with the
basis of the
P1, the additional basis of the
P2 then the additionnal basis of the
P6 (not with the basis of the
P1,
the additional basis of the P3 then the
additionnal basis of the P6, this is possible to
build the latter with "FEM_GEN_HIERARCHICAL(a,b))
The principal interest of the composite elements is to build
hierarchical elements. But this tool can also be used to build
piecewise polynomial elements.
composite element
"FEM_STRUCTURED_COMPOSITE(FEM_PK(2,1), 3)"
| Composition of a finite element
method on a element with S subdivisions |
| "FEM_STRUCTURED_COMPOSITE(FEM1,
S)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| degree of FEM1 |
dimension of FEM1 |
variable |
variable |
No (Q = 1) |
If FEM1 is |
piecewise |
It is important to use a corresponding composite integration
method.
hierarchical composite element
"FEM_PK_HIERARCHICAL_COMPOSITE(2,1,3)"
| Hierarchical composition of a
PK finite element method on a simplex with
S subdivisions |
|
"FEM_PK_HIERARCHICAL_COMPOSITE(P,K,S)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K |
P |
((SK+P)!)/((SK)! P!) |
variable |
No (Q = 1) |
If FEM1 is |
piecewise |
| hierarchical composition of a
hierarchical PK finite element method on a
simplex with S subdivisions |
|
"FEM_PK_FULL_HIERARCHICAL_COMPOSITE(P,K,S)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K |
P |
((SK+P)!)/((SK)! P!) |
variable |
No (Q = 1) |
If FEM1 is |
piecewise |
Other constructions are possible thanks to
"FEM_GEN_HIERARCHICAL(FEM1, FEM2)" and
"FEM_STRUCTURED_COMPOSITE(FEM1, S)"
It is important to use a corresponding composite integration
method.
RT0 elements in dimension two and three. (P+1 dof,
H(div))
| Raviart-Thomas of lowest order
element on simplices |
| "FEM_RT0(P)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 1 |
P |
P+1 |
H(div) |
Yes (Q = P) |
No |
Yes |
| Raviart-Thomas of lowest order
element on parallelepipeds (quadrilaterals, hexahedrals) |
| "FEM_RT0Q(P)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 1 |
P |
2P |
H(div) |
Yes (Q = P) |
No |
Yes |
Nedelec edge elements in dimension two and three.
(P(P+1)/2 dof, H(rot))
| Nedelec (or Whitney) edge
element |
| "FEM_NEDELEC(P)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 1 |
P |
P(P+1)/2 |
H(rot) |
Yes (Q = P) |
No |
Yes |
The 1D GaussLobatto PK element is similar to
the classical PK fem on the segment, but the
nodes are given by the Gauss-Lobatto-Legendre quadrature rule of
order 2K-1. This FEM is known to lead to better conditioned
linear systems, and can be used with the correspounding quadrature
to perform mass-lumping (on segments or parallelepipeds).
The polynomials coefficients have been pre-computed with Maple
(they require the inversion of an ill-conditionned system), hence
they are only available for the following values
of K: 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 16,
24, 32. Note that for K=1 and K=2, this is the
classical P1 and P2 fem.
| GaussLobatto PK
element on the segment |
|
"FEM_PK_GAUSSLOBATTO1D(K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| K |
1 |
K+1 |
C0 |
No (Q = 1) |
Yes |
Yes |
P3 Hermite element on a segment, 4
d.o.f., C1
Base functions on the reference element
| φ'0 =
(2x+1)(x-1)2, |
φ'1 =
x(x-1)2, |
| φ'2 =
x2(3-2x), |
φ'3 = x2(x -
1). |
This element is close to be
τ-equivalent but it is not.
On the real element the value of the gradient on vertices will be
multiplied by the gradient of the geometric transformation. The
matrix
M is not equal to identity but is still
diagonal.
| Hermite element on the
segment |
| "FEM_HERMITE(1)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
1 |
4 |
C1 |
No (Q = 1) |
No |
Yes |
P1 Lagrange element on a segment with
additional internal bubble function, 3 d.o.f., C0
| Lagrange P1
element with an additional internal bubble function |
| "FEM_PK_WITH_CUBIC_BUBBLE(1,
1)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 2 |
1 |
3 |
C0 |
No (Q = 1) |
Yes |
Yes |
 |
 |
| P1 with additional
bubble function, 4 d.o.f., C0 |
P2 with additional
bubble function, 7 d.o.f., C0 |
Lagrange element on a triangle with additional
internal bubble function
| Lagrange P1 or
P2 element with an additional internal bubble
function |
| "FEM_PK_WITH_CUBIC_BUBBLE(2,
K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
4 or 7 |
C0 |
No (Q = 1) |
Yes |
Yes |

P1 with additional bubble function, 4 d.o.f.,
C0
P1 Lagrange element on a triangle with
additional internal piecewise linear bubble function
| Lagrange P1
with an additional internal piecewise linear bubble
function |
|
"FEM_P1_PIECEWISE_LINEAR_BUBBLE" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 1 |
2 |
4 |
C0 |
No (Q = 1) |
Yes |
Piecewise |
P1 Lagrange element on a triangle with
additional bubble function on face 0, 4 d.o.f., C0
| Lagrange P1
element with an additional bubble function on face 0 |
|
"FEM_P1_BUBBLE_FACE(2)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 2 |
2 |
4 |
C0 |
No (Q = 1) |
Yes |
Yes |
P1 Lagrange element on a triangle with
additional d.o.f on face 0, 4 d.o.f., C0
| P1 Lagrange
element on a triangle with additional d.o.f on face 0 |
|
"FEM_P1_BUBBLE_FACE_LAG" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 2 |
2 |
4 |
C0 |
No (Q = 1) |
Yes |
Yes |
P1 non-conforming element on a triangle,
3 d.o.f., discontinuous
| P1
non-conforming element on a triangle |
| "FEM_P1_NONCONFORMING" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 1 |
2 |
3 |
discontinuous |
No (Q = 1) |
Yes |
Yes |
Hermite element on a triangle,
P3, 10 d.o.f., C0
Base functions on the reference element:
| φ'0 =
(1-x-y)(1+x+y-2x2-2y2-11xy), |
(φ'0(0,0) = 1), |
| φ'1 =
x(1-x-y)(1-x-2y), |
(∂xφ'1(0,0) = 1), |
| φ'2 =
y(1-x-y)(1-2x-y), |
(∂yφ'2(0,0) = 1), |
| φ'3 = -2x3 +
7 x2y + 7xy2 + 3x2 - 7xy, |
(φ'3(1,0) = 1), |
| φ'4 =
x3-2x2y-2xy2-x2+2xy, |
(∂xφ'4(1,0) = 1), |
| φ'5 = xy(y+2x-1), |
(∂yφ'5(1,0) = 1), |
| φ'6 = 7x2y +
7xy2 - 2y3+3y2-7xy, |
(φ'6(0,1) = 1), |
| φ'7 = xy(x+2y-1), |
(∂xφ'7(0,1) = 1), |
| φ'8 =
y3-2x2y-2xy2-y2+2xy, |
(∂yφ'8(0,1) = 1), |
| φ'9 = 27xy(1-x-y), |
(φ'9(1/3,1/3) =
1), |
|
This element is not
τ-equivalent (The matrix
M is
not equal to identity). On the real element linear combinaisons of
φ'4 and
φ'7 are used to
match the gradient on the corresponding vertex. Idem for the two
couples (
φ'5,
φ'8) and
(
φ'6,
φ'9) for the two
other vertices.
| Hermite element on a
triangle |
| "FEM_HERMITE(2)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
10 |
C0 |
No (Q = 1) |
No |
Yes |
triangle Morley element, P2, 6
d.o.f., C0
This element is not τ-equivalent (The matrix M
is not equal to identity). In particular, it can be used for
non-conforming discretization of fourth order problems, despite the
fact that it is not C0.
| Morley element on a
triangle |
| "FEM_MORLEY" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 2 |
2 |
6 |
|
No (Q = 1) |
No |
Yes |
Argyris element, P5, 21 d.o.f.,
C1
The base functions on the reference element are:
| φ'0(x,y) = 1 -
10x3 - 10y3 + 15x4 -
30x2y2 + 15y4 - 6x5 +
30x3y2 + 30x2y3 -
6y5, |
(φ'0(0,0) = 1), |
| φ'1(x,y) = x -
6x3 - 11xy2 + 8x4 +
10x2y2 + 18xy3 - 3x5 +
x3y2 - 10x2y3 -
8xy4, |
(∂xφ'1(0,0) = 1), |
| φ'2(x,y) = y -
11x2y - 6y3 + 18x3y +
10x2y2 + 8y4 - 8x4y -
10x3y2 + x2y3 -
3y5, |
(∂yφ'2(0,0) = 1), |
| φ'3(x,y) =
0.5x2 - 1.5x3 + 1.5x4 -
1.5x2y2 - 0.5x5 +
1.5x3y2 + x2y3, |
(∂2xxφ'3(0,0) = 1), |
| φ'4(x,y) = xy -
4x2y - 4xy2 + 5x3y +
10x2y2 + 5xy3 - 2x4y -
6x3y2 - 6x2y3 -
2xy4, |
(∂2xyφ'4(0,0) = 1), |
| φ'5(x,y) =
0.5y2 - 1.5y3 - 1.5x2y2
+ 1.5y4 + x3y2 +
1.5x2y3 - 0.5y5, |
(∂2yyφ'5(0,0) = 1), |
| φ'6(x,y) =
10x3 - 15x4 + 15x2y2 +
6x5 - 15x3y2 -
15x2y3, |
(φ'6(1,0) = 1), |
| φ'7(x,y) =
-4x3 + 7x4 - 3.5x2y2 -
3x5 + 3.5x3y2 +
3.5x2y3, |
(∂xφ'7(1,0) = 1), |
| φ'8(x,y) =
-5x2y + 14x3y +
18.5x2y2 - 8x4y -
18.5x3y2 -
13.5x2y3, |
(∂yφ'8(1,0) = 1), |
| φ'9(x,y) =
0.5x3 - x4 + 0.25x2y2 +
0.5x5 - 0.25x3y2 -
0.25x2y3, |
(∂2xxφ'9(1,0) = 1), |
| φ'10(x,y) =
x2y - 3x3y - 3.5x2y2 +
2x4y + 3.5x3y2 +
2.5x2y3, |
(∂2xyφ'10(1,0) = 1), |
| φ'11(x,y) =
1.25x2y2 - 0.75x3y2 -
1.25x2y3, |
(∂2yyφ'11(1,0) = 1), |
| φ'12(x,y) =
10y3 + 15x2y2 - 15y4 -
15x3y2 - 15x2y3 +
6y5, |
(φ'12(0,1) = 1), |
| φ'13(x,y) =
-5xy2 + 18.5x2y2 +
14xy3 - 13.5x3y2 -
18.5x2y3 - 8xy4, |
(∂xφ'13(0,1) = 1), |
| φ'14(x,y) =
-4y3 - 3.5x2y2 + 7y4 +
3.5x3y2 + 3.5x2y3 -
3y5, |
(∂yφ'14(0,0) = 1), |
| φ'15(x,y) =
1.25x2y2 - 1.25x3y2 -
0.75x2y3, |
(∂2xxφ'15(0,1) = 1), |
| φ'16(x,y) =
xy2 - 3.5x2y2 - 3xy3 +
2.5x3y2 + 3.5x2y3 +
2xy4, |
(∂2xyφ'16(0,1) = 1), |
| φ'17(x,y) =
0.5y3 + 0.25x2y2 - y4 -
0.25x3y2 - 0.25x2y3 +
0.5y5, |
(∂2yyφ'17(0,1) = 1), |
| φ'18(x,y) =
sqrt(2)(-8x2y2 + 8x3y2
+ 8x2y3), |
(sqrt(0.5)(∂xφ'18(0.5,0.5) +
∂yφ'18(0.5,0.5)) = 1), |
| φ'19(x,y) =
-16xy2 + 32x2y2 + 32xy3
- 16x3y2 - 32x2y3 -
16xy4, |
(-∂xφ'19(0,0.5) = 1), |
| φ'20(x,y) =
-16x2y + 32x3y + 32x2y2
- 16x4y - 32x3y2 -
16x2y3, |
(-∂yφ'20(0.5,0) = 1), |
|
This element is not
τ-equivalent (The matrix
M is
not equal to identity). On the real element linear combinaisons of
the transformed base functions
φ'i are used
to match the gradient, the second derivatives and the normal
derivatives on the faces. Note that the use of the matrix
M
allows to define Argyris element even with nonlinear geometric
transformations (for instance to treat curved
boundaries).
| Argyris element on a
triangle |
| "FEM_ARGYRIS" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 5 |
2 |
21 |
C1 |
No (Q = 1) |
No |
Yes |
Hsieh-Clough-Tocher (HCT) element,
P3, 12 d.o.f., C1
This element is not τ-equivalent. This is a composite
element. Polynomial of degree 3 on each of the three sub-triangles
(see figure * and [1]). It is strongly advised
to use a IM_HCT_COMPOSITE
integration method with this finite element. The numeration of the
dof is the following : 0, 3 and 6 for the lagrange dof on the first
second and third vertex respectively; 1, 4, 7 for the derivative
with respects to the first variable; 2, 5, 8 for the derivative
with respects to the second variable and 9, 10, 11 for the normal
derivatives on face 0, 1, 2 respectively.
| HCT element on a
triangle |
| "FEM_HCT_TRIANGLE" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
12 |
C1 |
No (Q = 1) |
No |
piecewise |
Reduced Hsieh-Clough-Tocher (reduced HCT)
element, P3, 9 d.o.f., C1
This element exists also in its reduced form, where the normal
derivatives is assumed to be polynomial of degree one on each edge
(see figure *)
| Reduced HCT element on a
triangle |
|
"FEM_REDUCED_HCT_TRIANGLE" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
9 |
C1 |
No (Q = 1) |
No |
piecewise |
Composite element on quadrilaterals,
piecewise P3, 16 d.o.f., C1
This element is not τ-equivalent. This is a composite
element. Polynomial of degree 3 on each of the four sub-triangles
(see figure *). At least on the reference
element it correponds to the Fraeijs de Veubeke-Sander element (see
[1]). It is strongly
advised to use a IM_QUADC1_COMPOSITE
integration method with this finite element.
| C1 composite
element on a quadrilateral (FVS) |
| "FEM_QUADC1_COMPOSITE" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
16 |
C1 |
No (Q = 1) |
No |
piecewise |
Reduced composite element on quadrilaterals,
piecewise P3, 12 d.o.f., C1
This element exists also in its reduced form, where the normal
derivatives is assumed to be polynomial of degree one on each edge
(see figure *)
| Reduced C1
composite element on a quadrilateral (reduced FVS) |
|
"FEM_REDUCED_QUADC1_COMPOSITE" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
2 |
12 |
C1 |
No (Q = 1) |
No |
piecewise |
 |
 |
 |
| P1 with additional
bubble function, 5 d.o.f., C0 |
P2 with additional
bubble function, 11 d.o.f., C0 |
P3 with additional
bubble function, 21 d.o.f., C0 |
Lagrange element on a tetrahedron with additional
internal bubble function.
| PK Lagrange
element with an additional internal bubble function |
| "FEM_PK_WITH_CUBIC_BUBBLE(3,
K)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 4 |
3 |
5, 11 or 21 |
C0 |
No (Q = 1) |
Yes |
Yes |
P1 Lagrange element on a tetrahedron
with additional bubble function on face 0, 5 d.o.f.,
C0
| Lagrange P1
element with an additional bubble function on face 0 |
|
"FEM_P1_BUBBLE_FACE(3)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
3 |
5 |
C0 |
No (Q = 1) |
Yes |
Yes |
Hermite element on a tetrahedron,
P3, 20 d.o.f., C0
Base functions on the reference element:
| φ'0(x,y) = 1 -
3x2 - 13xy - 13xz - 3y2 - 13yz -
3z2 + 2x3 + 13x2y +
13x2z |
|
|
+ 13xy2 + 33xyz + 13xz2 + 2y3 +
13y2z + 13yz2 + 2z3, |
(φ'0(0,0,0) = 1), |
| φ'1(x,y) = x -
2x2 - 3xy - 3xz + x3 + 3x2y +
3x2z + 2xy2 + 4xyz + 2xz2, |
(∂xφ'1(0,0,0) = 1), |
| φ'2(x,y) = y - 3xy -
2y2 - 3yz + 2x2y + 3xy2 + 4xyz +
y3 + 3y2z + 2yz2, |
(∂yφ'2(0,0,0) = 1), |
| φ'3(x,y) = z - 3xz -
3yz - 2z2 + 2x2z + 4xyz + 3xz2 +
2y2z + 3yz2 + z3, |
(∂zφ'3(0,0,0) = 1), |
| φ'4(x,y) =
3x2 - 7xy - 7xz - 2x3 + 7x2y +
7x2z + 7xy2 + 7xyz + 7xz2, |
(φ'4(1,0,0) = 1), |
| φ'5(x,y) =
-x2 + 2xy + 2xz + x3 - 2x2y -
2x2z - 2xy2 - 2xyz - 2xz2, |
(∂xφ'5(1,0,0) = 1), |
| φ'6(x,y) = -xy +
2x2y + xy2, |
(∂yφ'6(1,0,0) = 1), |
| φ'7(x,y) = -xz +
2x2z + xz2, |
(∂zφ'7(1,0,0) = 1), |
| φ'8(x,y) = -7xy +
3y2 - 7yz + 7x2y + 7xy2 + 7xyz -
2y3 + 7y2z + 7yz2, |
(φ'8(0,1,0) = 1), |
| φ'9(x,y) = -xy +
x2y + 2xy2, |
(∂xφ'9(0,1,0) = 1), |
| φ'10(x,y) = 2xy -
y2 + 2yz - 2x2y - 2xy2 - 2xyz +
y3 - 2y2z - 2yz2, |
(∂yφ'10(0,1,0) = 1), |
| φ'11(x,y) = -yz +
2y2z + yz2, |
(∂zφ'11(0,1,0) = 1), |
| φ'12(x,y) = -7xz - 7yz
+ 3z2 + 7x2z + 7xyz + 7xz2 +
7y2z + 7yz2 - 2z3, |
(φ'12(0,0,1) = 1), |
| φ'13(x,y) = -xz +
x2z + 2xz2, |
(∂xφ'13(0,0,1) = 1), |
| φ'14(x,y) = -yz +
y2z + 2yz2, |
(∂yφ'14(0,0,1) = 1), |
| φ'15(x,y) = 2xz + 2yz -
z2 - 2x2z - 2xyz - 2xz2 -
2y2z - 2yz2 + z3, |
(∂zφ'15(0,0,1) = 1), |
| φ'16(x,y) = 27xyz, |
(φ'16(1/3,1/3,1/3) =
1), |
| φ'17(x,y) = 27yz -
27xyz - 27y2z - 27yz2, |
(φ'17(0,1/3,1/3) =
1), |
| φ'18(x,y) = 27xz -
27x2z - 27xyz - 27xz2, |
(φ'18(1/3,0,1/3) =
1), |
| φ'19(x,y) = 27xy -
27x2y - 27xy2 - 27xyz, |
(φ'19(1/3,1/3,0) =
1), |
|
This element is not
τ-equivalent (The matrix
M is
not equal to identity). On the real element linear combinaisons of
φ'8,
φ'12 and
φ'16 are used to match the gradient on the
corresponding vertex. Idem on the orther vertices.
| Hermite element on a
tetrahedron |
| "FEM_HERMITE(3)" |
| Degree |
dimension |
d.o.f. number |
class |
vectorial |
τ-equivalent |
Polynomial |
| 3 |
3 |
20 |
C0 |
No (Q = 1) |
No |
Yes |