Manpages

NAME

integrate_option - expression integration options (rheolef-7.2)

DESCRIPTION

This class sends options to the integrate(3) function, when building a form(2) or a field(2).

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(1) 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.

QUADRATURE FAMILIES

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(2). 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 order parameter refers to the number of nodes used and not to the degree of polynomial exactly integated. For instance, order=1 refers to the trapezoidal formulae and for the general case, there are order+1 nodes per edges.

BOOLEAN FLAGS

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 r weight in the r dr dz measure during the numerical integration performed the integrate(3) function. This feature is useful for computing the stream function in the axisymmetric case. Default is ignore_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 P1 element. Default is lump=false.

IMPLEMENTATION

This documentation has been generated from file fem/geo_element/integrate_option.h

class integrate_option {
public:
// typedefs:

typedef size_t size_type;

typedef enum {
gauss = 0,
gauss_lobatto = 1,
gauss_radau = 2,
middle_edge = 3,
superconvergent = 4,
equispaced = 5,
max_family = 6
} family_type; // update also family_name[] in quatrature.cc

static const size_type unset_order = std::numeric_limits<size_type>::max();
static const size_type default_order = unset_order;
static const family_type default_family = gauss;

// allocators:

integrate_option(
family_type ft = default_family,
size_type k = default_order);

integrate_option (const std::string& name);
integrate_option (const integrate_option& iopt);
integrate_option& operator= (const integrate_option& iopt);

// accessors & modifiers:

std::string name() const;
size_t get_order() const;
family_type get_family() const;
std::string get_family_name() const;
void reset (const std::string& name);
void set_order (size_t r);
void set_family (family_type type);
void set_family (std::string name);

// data:

bool invert, ignore_sys_coord, lump;

};

AUTHOR

Pierre Saramito <Pierre.Saramito [AT] imag.fr>

COPYRIGHT

Copyright (C) 2000-2018 Pierre Saramito <Pierre.Saramito [AT] imag.fr> GPLv3+: GNU GPL version 3 or later <http://gnu.org/licenses/gpl.html>;. This is free software: you are free to change and redistribute it. There is NO WARRANTY, to the extent permitted by law.