.TH "integrate_option" 3rheolef "Version 7.2" "rheolef" \" -*- nroff -*- .ad l .nh .SH NAME integrate_option \- expression integration options (rheolef-7\&.2) .PP .SH "DESCRIPTION" .PP This class sends options to the \fBintegrate(3)\fP function, when building a \fBform(2)\fP or a \fBfield(2)\fP\&. .PP 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 \fCorder\fP of the quadrature formula\&. See also \fBquadrature(1)\fP for examples of quadrature formula\&. .PP In addition to the customization of the quadrature formula, the present class provides some booleaan flags, e\&.g\&. computation of the inverse matrix\&. .SH "QUADRATURE FAMILIES" .PP \fCgauss\fP .PP .RS 4 The Gauss formula is the default\&. .RE .PP \fCgauss_lobatto\fP .PP .RS 4 The Gauss-Lobatto formula is useful in some cases, e\&.g\&. when using the method of \fBcharacteristic(2)\fP\&. It is actually implemented for order <= 2\&. .RE .PP \fCgauss_radau\fP .PP .RS 4 The Gauss-Radau formula is another classical alternative\&. .RE .PP \fCmiddle_edge\fP .PP .RS 4 The nodes are located at the midle edge\&. This formula is useful in some special cases and its order is limited\&. .RE .PP \fCsuperconvergent\fP .PP .RS 4 Another special case, for testing superconvergence of the finite element method at some specific nodes\&. .RE .PP \fCequispaced\fP .PP .RS 4 When the solution has low regularity, e\&.g\&. is discontinuous across elements, this simple formula is also useful\&. In that case, the \fCorder\fP parameter refers to the number of nodes used and not to the degree of polynomial exactly integated\&. For instance, \fCorder=1\fP refers to the trapezoidal formulae and for the general case, there are \fCorder+1\fP nodes per edges\&. .RE .PP .SH "BOOLEAN FLAGS" .PP \fCinvert\fP .PP .RS 4 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 \fCinvert=false\fP\&. .RE .PP \fCignore_sys_coord\fP .PP .RS 4 This flag has effects only for axisymmetric coordinate systems\&. When set, it omits the \fCr\fP weight in the \fCr dr dz\fP measure during the numerical integration performed the \fBintegrate(3)\fP function\&. This feature is useful for computing the stream function in the axisymmetric case\&. Default is \fCignore_sys_coord=false\fP\&. .RE .PP \fClump\fP .PP .RS 4 This flag, when set, performs a \fImass lumping\fP procedure on the matrix at the element level: .RE .PP a(i,i) += sum(j!=i) a(i,j) .PP .RS 4 The resulting matrix is diagonal\&. This feature is useful for computing a diagonal approximation of the mass matrix for the continuous \fCP1\fP element\&. Default is \fClump=false\fP\&. .RE .PP .SH "IMPLEMENTATION" .PP This documentation has been generated from file fem/geo_element/integrate_option\&.h .PP .PP .nf 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::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; .fi .PP .PP .nf }; .fi .PP .SH AUTHOR Pierre Saramito .SH COPYRIGHT Copyright (C) 2000-2018 Pierre Saramito GPLv3+: GNU GPL version 3 or later . This is free software: you are free to change and redistribute it. There is NO WARRANTY, to the extent permitted by law.