|
Rheolef
7.2
an efficient C++ finite element environment
|
expression integration options
This class sends options to the integrate function, when building a form or a field.
Its allows one to choose the quadrature formula used during the numerical integration (Gauss, Gauss-Lobatto, etc) and the polynomial degree that is exactly integrated. This exactly integrated polynomial degree is called here the order of the quadrature formula. See also quadrature for examples of quadrature formula.
In addition to the customization of the quadrature formula, the present class provides some booleaan flags, e.g. computation of the inverse matrix.
gauss
The Gauss formula is the default.
gauss_lobatto
The Gauss-Lobatto formula is useful in some cases, e.g. when using the method of
characteristic. It is actually implemented for order <= 2.
gauss_radau
The Gauss-Radau formula is another classical alternative.
middle_edge
The nodes are located at the midle edge. This formula is useful in some special cases and its order is limited.
superconvergent
Another special case, for testing superconvergence of the finite element method at some specific nodes.
equispaced
When the solution has low regularity, e.g. is discontinuous across elements, this simple formula is also useful. In that case, the
orderparameter refers to the number of nodes used and not to the degree of polynomial exactly integated. For instance,order=1refers to the trapezoidal formulae and for the general case, there areorder+1nodes per edges.
invert
This flag, when set, performs a local inversion on the matrix at the element level: This procedure is allowed only when the global matrix is block diagonal, e.g. for discontinuous or bubble approximations. This property is true when basis functions have a compact support inside exactly one element. Default is
invert=false.
ignore_sys_coord
This flag has effects only for axisymmetric coordinate systems. When set, it omits the
rweight in ther dr dzmeasure during the numerical integration performed theintegratefunction. This feature is useful for computing the stream function in the axisymmetric case. Default isignore_sys_coord=false.
lump
This flag, when set, performs a mass lumping procedure on the matrix at the element level:
a(i,i) += sum(j!=i) a(i,j)
The resulting matrix is diagonal. This feature is useful for computing a diagonal approximation of the mass matrix for the continuous
P1element. Default islump=false.
This documentation has been generated from file fem/geo_element/integrate_option.h