Appendix B. Cubature method list

The integration methods are of two kinds. Exact integrations of polynomials and approximated integrations (cubature formulas) of any function. The exact integration can only be used if all the elements are polynomial and if the geometric transformation is linear.

A descriptor on an integration method is given by the function

ppi = getfem::int_method_descriptor("name of method");

where "name of method" is a string to be choosen among the existing methods.

The program integration located in the tests directory lists and checks the degree of each integration method.

Exact Integration methods

The list of available Exact integration methods is the following

"IM_NONE()" Dummy integration method.
"IM_EXACT_SIMPLEX(n)" Description of the exact integration of polynomials on the simplex of reference of dimension n.
"IM_PRODUCT(a, b)" Description of the exact integration on the convex which is the direct product of the convex in a and in b.
"IM_EXACT_PARALLELEPIPED(n)" Description of the exact integration of polynomials on the parallelepiped of reference of dimension n
"IM_EXACT_PRISM(n)" Description of the exact integration of polynomials on the prism of reference of dimension n
Even though a description of exact integration method exists on parallelepipeds or prisms, most of the time the geometric transformations on such elements are not linear and the exact integration cannot be used.
Beware : In fact a lot of computation cannot be done with exact integration methods. So, it is recommended to use cubature formulas instead.

Newton cotes Integration methods

use "IM_NC(N,K)", "IM_NC_PARALLELEPIPED(N,K)" and "IM_NC_PRISM(N,K)" to have the Newton cotes integration of order K on simplices, parallelepipeds and prisms respectively.

Gauss Integration methods on dimension 1

use "IM_GAUSS1D(K)" to have the Gauss-Legendre integration on the segment of order K (with K/2 + 1 points), and "IM_GAUSSLOBATTO1D(K)" to have the Gauss-Lobatto-Legendre integration on the segment of order K (with K/2 + 1 points). The latter integration method is only available for odd values of K. The Gauss-Lobatto integration method can be used in conjunction with "FEM_PK_GAUSSLOBATTO1D(K/2)" to perform mass-lumping.

Gauss Integration methods on dimension 2

graphic coordinates
x y
weights function to call / order
IM_TRIANGLE(1)
1/3 1/3
1/2
"IM_TRIANGLE(1)" 1 point, order 1.
1/6 1/6
2/3 1/6
1/6 2/3
1/6
1/6
1/6
"IM_TRIANGLE(2)" 3 points, order 2.
1/3 1/3
1/5 1/5
3/5 1/5
1/5 3/5
-27/96
25/96
25/96
25/96
"IM_TRIANGLE(3)" 4 points, order 3.
a a
1-2a a
a 1-2a
b b
1-2b b
b 1-2b
c
c
c
d
d
d
"IM_TRIANGLE(4)" 6 points, order 4, a = 0.445948490915965, b = 0.091576213509771, c = 0.111690794839005, d = 0.054975871827661.
1/3 1/3
a a
1-2a a
a 1-2a
b b
1-2b b
b 1-2b
9/80
c
c
c
d
d
d
"IM_TRIANGLE(5)" 7 points, order 5, a = (6+sqrt(15))/(21), b = 4/7 - a, c = (155+sqrt(15))/(2400), d = 31/240 - c.
a a
1-2a a
a 1-2a
b b
1-2b b
b 1-2b
c d
d c
1-c-d c
1-c-d d
c 1-c-d
d 1-c-d
e
e
e
f
f
f
g
g
g
g
g
g
"IM_TRIANGLE(6)" 12 points, order 6, a = 0.063089104491502, b = 0.249286745170910, c = 0.310352451033785, d = 0.053145049844816, e = 0.025422453185103, f = 0.058393137863189, g = 0.041425537809187.
a a
b a
a b
c e
d c
e d
d e
c d
e c
f f
g f
f g
1/3 1/3
h
h
h
i
i
i
i
i
i
j
j
j
k
"IM_TRIANGLE(7)" 13 points, order 7, a = 0.0651301029022, b = 0.8697397941956, c = 0.3128654960049, d = 0.6384441885698, e = 0.0486903154253, f = 0.2603459660790, g = 0.4793080678419, h = 0.0266736178044, i = 0.0385568804451, j = 0.0878076287166, k = -0.0747850222338.
"IM_TRIANGLE(8)" (see [3]) 16 points, order 8
"IM_TRIANGLE(9)" (see [3]) 19 points, order 9
"IM_TRIANGLE(10)" (see [3]) 25 points, order 10
"IM_TRIANGLE(13)" (see [3]) 37 points, order 13
1/2+sqrt(1/6) 1/2
1/2-sqrt(1/24) 1/2±sqrt(1/8)
1/3
1/3
"IM_QUAD(2)" 3 points, order 2.
1/2±sqrt(1/6) 1/2
1/2 1/2±sqrt(1/6)
1/4
1/4
"IM_QUAD(3)" 4 points, order 3.
1/2 1/2
1/2 ±sqrt(7/30) 1/2
1/2±sqrt(1/12) 1/2±sqrt(3/20)
2/7
5/63
5/36
"IM_QUAD(5)" 7 points, order 5.
"IM_QUAD(7)" 12 points, order 7
"IM_QUAD(9)" 20 points, order 9
"IM_QUAD(17)" 70 points, order 17
 

There is also the IM_GAUSS_PARALLELEPIPED(n,k) which is a direct product of 1D gauss integrations.
Important note: do not forget that IM_QUAD(k) is exact for polynomials up to degree k, and that a Qk polynomial has a degree of 2*k. For example, IM_QUAD(7) cannot integrate exactly the product of two Q2 polynomials. On the other hand, IM_GAUSS_PARALLELEPIPED(2,4) can integrate exactly that product...

Gauss Integration methods on dimension 3

graphic coordinates
x y z
weights function to call / order
1/4 1/4 1/4
1/6
"IM_TETRAHEDRON(1)" 1 point, order 1.
a a a
a b a
a a b
b a a
1/24
1/24
1/24
1/24
"IM_TETRAHEDRON(2)" 4 points, order 2 a = (5 - sqrt(5))/(20), b = (5 + 3sqrt(5))/(20).
1/4 1/4 1/4
1/6 1/6 1/6
1/6 1/2 1/6
1/6 1/6 1/2
1/2 1/6 1/6
-2/15
3/40
3/40
3/40
"IM_TETRAHEDRON(3)" 5 points, order 3
1/4 1/4 1/4
a a a
a a c
a c a
c a a
b b b
b b d
b d b
d b b
e e f
e f e
f e e
e f f
f e f
f f e
8/405
h
h
h
h
i
i
i
i
5/567
5/567
5/567
5/567
5/567
5/567
"IM_TETRAHEDRON(5)" 15 points, order 5 a = (7 + sqrt(15))/(34), b = (7 - sqrt(15))/(34), c = (13 + 3sqrt(15))/(34), d = (13 - 3sqrt(15))/(34), e = (5 - sqrt(15))/(20), f = (5 + sqrt(15))/(20), h = (2665 - 14sqrt(15))/(226800), i = (2665 + 14sqrt(15))/(226800),

Others methods are:

name convex type nb of points
IM_TETRAHEDRON(6) 3D simplex 24
IM_TETRAHEDRON(8) 3D simplex 43
IM_SIMPLEX4D(3) 4D simplex 6
IM_HEXAHEDRON(5) 3D parallelepipeded 14
IM_HEXAHEDRON(9) 3D parallelepipeded 58
IM_HEXAHEDRON(11) 3D parallelepipeded 90
IM_CUBE4D(5) 4D parallelepipeded 24
IM_CUBE4D(9) 4D parallelepipeded 145

Direct product of integration methods

You can use "IM_PRODUCT(IM1, IM2)" to produce integration methods on quadrilateral or prisms. It gives the direct product of two integration mathods. For instance IM_GAUSS_PARALLELEPIPED(2,k) is an alias for IM_PRODUCT(IM_GAUSS1D(2,k),IM_GAUSS1D(2,k)) and ca be use instead of the IM_QUAD integrations.

Composite integration methods

composite method "IM_STRUCTURED_COMPOSITE(IM_TRIANGLE(2), 3)"
 

use "IM_STRUCTURED_COMPOSITE(IM1, S)" to copy IM1 on an element with S subdivisions. The resulting integration method has the same order but with more points. This could be more stable to use composite method rather than to improve the order of the method. Those methods have to be used also with composite elements. Most of the time for composite element, it is preferable to choose the basic method IM1 with no points on the boundary (because the gradient coulb be not defined on the boundary of sub-elements).

For the HCT element, it is advised to use the IM_HCT_COMPOSITE(im) composite integration (which split the original triangle into 3 sub-triangles).