.TH "solver" 4rheolef "Sat Mar 13 2021" "Version 7.1" "rheolef" \" -*- nroff -*- .ad l .nh .SH NAME solver \- direct and iterative solver (rheolef-7\&.1) .PP .SH "DESCRIPTION" .PP The class implements the numerical resolution of a linear system\&. Let \fCa\fP be a square and invertible matrix in the \fBcsr(4)\fP sparse format\&. The construction of a solver writes: .PP .nf solver sa (a); .fi .PP and the resolution of \fCa*x = b\fP expresses simply: .PP .nf vec x = sa.solve(b); .fi .PP When the matrix is modified in a computation loop, the solver could be re-initialized by: .PP .nf sa.update_values (new_a); vec x = sa.solve(b); .fi .PP .SH "DIRECT VERSUS ITERATIVE" .PP The choice between a direct or an iterative method for solving the linear system is by default performed automatically: it depends upon the sparsity pattern of the matrix, in order to achieve the best performances\&. The \fBsolver_option(4)\fP class allows one to change this default behavior\&. .PP .nf solver_option sopt; sopt.iterative = true; solver sa (a, sopt); vec x = sa.solve(b); .fi .PP The direct approach bases on the Choleski factorization for a symmetric definite positive matrix, and on the LU one otherwise\&. Conversely, the iterative approach bases on the \fBcg(5)\fP conjugate gradient algorithm for a symmetric definite positive matrix, and on the \fBgmres(5)\fP algorithm otherwise\&. .SH "COMPUTING THE DETERMINANT" .PP This feature is useful e\&.g\&. when tracking a change of sign in the determinant of a matrix, e\&.g\&. during the \fBcontinuation(3)\fP algorithm\&. When using a direct method, the determinant of the matrix can be computed as: .PP .nf solver_option sopt; sopt.iterative = false; solver sa (a, sopt); cout << sa.det().mantissa << "*" << sa.det().base << "^" << sa.det().exponant << endl; .fi .PP The \fCsa\&.det()\fP method returns an object of type \fCsolver::determinant_type\fP that contains a mantissa and an exponent in a given base (generally 2 or 10)\&. For some rare direct solvers, the computation of the determinant is not yet fully supported: it is the case for the Cholesky factorization from the \fCeigen\fP library\&. In you find such a problem, please switch to another solver library, see the \fBsolver_option(4)\fP class\&. .SH "AUTOMATIC CHOICE AND CUSTOMIZATION" .PP When the matrix is obtained from the finite element discretization of 3D partial differential problems, the iterative solver is the default choice\&. Otherwise, the direct solver is selected\&. .PP More precisely, the choice between direct or iterative solver depends upon the \fCa\&.pattern_dimension()\fP property of the \fBcsr(4)\fP sparse matrix\&. When this pattern dimension is 3, an iterative method is faster and less memory consuming than a direct one\&. See \fBusersguide\fP for a discussion on this subject\&. .PP The symmetry-positive-definiteness of the matrix is tested via the \fCa\&.is_symmetric()\fP and \fCa\&.is_definite_positive()\fP properties of the \fBcsr(4)\fP sparse matrix\&. These properties determine the choices between Cholesky/LU methods for the direct case, and between \fCcg/gmres\fP algorithms for the iterative one\&. Most of the time, these properties are automatically well initialized by the finite element assembly procedure, via the \fBintegrate(3)\fP function\&. .PP Nevertheless, in some special cases, e\&.g\&. a linear combination of matrix, or when the matrix has been read from a file, it could be necessary to force either the symmetry or the positive-definiteness property by the appropriate \fBcsr(4)\fP member function before to send the matrix to the solver\&. .SH "PRECONDITIONNERS FOR ITERATIVE SOLVERS" .PP When using an iterative method, the default is to perform no preconditionning\&. Several preconditionners are available: the \fBmic(5)\fP modified incomplete Cholesky for symmetric matrix and the \fBilut(5)\fP incomplete LU one for unsymmetric matrix and the do-nothing \fBeye(5)\fP identity preconditionner\&. A preconditionner can be supplied via: .PP .nf solver_option sopt; sopt.iterative = true; solver sa (a, sopt); sa.set_preconditionner (ilut(a)); vec x = sa.solve(b); .fi .PP Note also the \fBeye(5)\fP that returns the solver for the identity matrix: it could be be used for specifying that we do not use a preconditionner\&. This is the default behavior\&. The \fCset_preconditionner\fP member function should be called before the first call to the \fCsolve\fP method: if no preconditioner has been defined at the first call to \fCsolve\fP, the default \fBeye(5)\fP preconditionner is selected\&. .SH "IMPLEMENTATION" .PP This documentation has been generated from file linalg/lib/solver\&.h .PP The \fCsolver\fP class is simply an alias to the \fC\fBsolver_basic\fP\fP class .PP .PP .nf typedef solver_basic solver; .fi .PP \fB\fP .RS 4 .RE .PP The \fC\fBsolver_basic\fP\fP class provides an interface to a data container: .PP .PP .nf template class solver_basic: public smart_pointer_clone > { public: // typedefs: typedef solver_abstract_rep rep; typedef smart_pointer_clone base; typedef typename rep::size_type size_type; typedef typename rep::determinant_type determinant_type; // allocators: solver_basic (); explicit solver_basic (const csr& a, const solver_option& opt = solver_option()); void update_values (const csr& a); // accessors: vec trans_solve (const vec& b) const; vec solve (const vec& b) const; determinant_type det() const; const solver_option& option() const; void set_preconditioner (const solver_basic&); bool initialized() const; std::string name() const; }; .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.