|
| MINRES () |
|
template<typename MatrixDerived > |
| MINRES (const EigenBase< MatrixDerived > &A) |
|
| ~MINRES () |
|
template<typename Rhs , typename Dest > |
void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
|
template<typename Rhs , typename Dest > |
void | _solve_impl (const Rhs &b, Dest &x) const |
|
| IterativeSolverBase () |
|
| IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
|
| IterativeSolverBase (IterativeSolverBase &&)=default |
|
| ~IterativeSolverBase () |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | factorize (const EigenBase< MatrixDerived > &A) |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | compute (const EigenBase< MatrixDerived > &A) |
|
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
RealScalar | tolerance () const |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | setTolerance (const RealScalar &tolerance) |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
Index | maxIterations () const |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | setMaxIterations (Index maxIters) |
|
Index | iterations () const |
|
RealScalar | error () const |
|
const SolveWithGuess< MINRES< MatrixType_, UpLo_, Preconditioner_ >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
ComputationInfo | info () const |
|
void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime !=1 &&DestDerived::ColsAtCompileTime !=1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
|
void | _solve_impl (const Rhs &b, Dest &x) const |
|
MINRES< MatrixType_, UpLo_, Preconditioner_ > & | derived () |
|
const MINRES< MatrixType_, UpLo_, Preconditioner_ > & | derived () const |
|
| SparseSolverBase () |
|
| SparseSolverBase (SparseSolverBase &&other) |
|
| ~SparseSolverBase () |
|
Derived & | derived () |
|
const Derived & | derived () const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
template<typename Rhs , typename Dest > |
void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
|
template<typename MatrixType_, int UpLo_, typename Preconditioner_>
class Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >
A minimal residual solver for sparse symmetric problems.
This class allows to solve for A.x = b sparse linear problems using the MINRES algorithm of Paige and Saunders (1975). The sparse matrix A must be symmetric (possibly indefinite). The vectors x and b can be either dense or sparse.
- Template Parameters
-
MatrixType_ | the type of the sparse matrix A, can be a dense or a sparse matrix. |
UpLo_ | the triangular part that will be used for the computations. It can be Lower, Upper, or Lower|Upper in which the full matrix entries will be considered. Default is Lower. |
Preconditioner_ | the type of the preconditioner. Default is DiagonalPreconditioner |
The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations and NumTraits<Scalar>::epsilon() for the tolerance.
This class can be used as the direct solver classes. Here is a typical usage example:
SparseMatrix<double>
A(
n,
n);
MINRES<SparseMatrix<double> > mr;
std::cout << "#iterations: " << mr.iterations() << std::endl;
std::cout << "estimated error: " << mr.error() << std::endl;
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
list x
Definition: plotDoE.py:28
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.
MINRES can also be used in a matrix-free context, see the following example .
- See also
- class ConjugateGradient, BiCGSTAB, SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
template<typename MatrixType_ , int UpLo_, typename Preconditioner_ >
template<typename Rhs , typename Dest >
void Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::_solve_vector_with_guess_impl |
( |
const Rhs & |
b, |
|
|
Dest & |
x |
|
) |
| const |
|
inline |
231 TransposeInput = (!MatrixWrapper::MatrixFree) && (
UpLo == (
Lower |
Upper)) && (!MatrixType::IsRowMajor) &&
234 typedef std::conditional_t<TransposeInput, Transpose<const ActualMatrixType>,
ActualMatrixType const&>
237 MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
238 typedef std::conditional_t<
UpLo == (
Lower |
Upper), RowMajorWrapper,
244 RowMajorWrapper row_mat(
matrix());
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
Definition: IterativeSolverBase.h:371
Index maxIterations() const
Definition: IterativeSolverBase.h:251
MatrixWrapper::ActualMatrixType ActualMatrixType
Definition: IterativeSolverBase.h:372
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:382
RealScalar m_tolerance
Definition: IterativeSolverBase.h:385
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
RealScalar m_error
Definition: IterativeSolverBase.h:387
Index m_iterations
Definition: IterativeSolverBase.h:388
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
constexpr bool check_implication(bool a, bool b)
Definition: Meta.h:740
EIGEN_DONT_INLINE void minres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
Definition: MINRES.h:32
Type
Type of JSON value.
Definition: rapidjson.h:513
@ IsComplex
Definition: NumTraits.h:176
References b, Eigen::internal::check_implication(), EIGEN_STATIC_ASSERT, Eigen::Lower, Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::m_error, Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::m_info, Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::m_iterations, Eigen::IterativeSolverBase< Derived >::m_preconditioner, Eigen::IterativeSolverBase< Derived >::m_tolerance, Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::matrix(), Eigen::IterativeSolverBase< Derived >::maxIterations(), Eigen::internal::minres(), Eigen::NoConvergence, Eigen::Success, Eigen::MINRES< MatrixType_, UpLo_, Preconditioner_ >::UpLo, Eigen::Upper, and plotDoE::x.