16#include "./InternalHeaderCheck.h"
31template <
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
32EIGEN_DONT_INLINE
void minres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond,
33 Index& iters,
typename Dest::RealScalar& tol_error) {
35 typedef typename Dest::RealScalar RealScalar;
36 typedef typename Dest::Scalar Scalar;
37 typedef Matrix<Scalar, Dynamic, 1> VectorType;
40 const RealScalar rhsNorm2(rhs.squaredNorm());
49 const Index maxIters(iters);
50 const Index N(mat.cols());
51 const RealScalar threshold2(tol_error * tol_error * rhsNorm2);
55 VectorType v(VectorType::Zero(N));
56 VectorType v_new(rhs - mat * x);
57 RealScalar residualNorm2(v_new.squaredNorm());
59 VectorType w_new(precond.solve(v_new));
61 RealScalar beta_new2(v_new.dot(w_new));
62 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
63 RealScalar beta_new(
sqrt(beta_new2));
64 const RealScalar beta_one(beta_new);
67 RealScalar c_old(1.0);
69 RealScalar s_old(0.0);
71 VectorType p_old(VectorType::Zero(N));
76 while (iters < maxIters) {
87 const RealScalar beta(beta_new);
93 v_new.noalias() = mat * w - beta * v_old;
94 const RealScalar alpha = v_new.dot(w);
96 w_new = precond.solve(v_new);
97 beta_new2 = v_new.dot(w_new);
98 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
99 beta_new =
sqrt(beta_new2);
102 const RealScalar r2 = s * alpha + c * c_old * beta;
103 const RealScalar r3 = s_old * beta;
104 const RealScalar r1_hat = c * alpha - c_old * s * beta;
105 const RealScalar r1 =
sqrt(std::pow(r1_hat, 2) + std::pow(beta_new, 2));
114 p.noalias() = (w - r2 * p_old - r3 * p_oold) / r1;
115 x += beta_one * c * eta * p;
119 residualNorm2 *= s * s;
121 if (residualNorm2 < threshold2) {
131 tol_error = std::sqrt(residualNorm2 / rhsNorm2);
136template <
typename MatrixType_,
int UpLo_ = Lower,
typename Preconditioner_ = IdentityPreconditioner>
141template <
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
142struct traits<MINRES<MatrixType_, UpLo_, Preconditioner_> > {
143 typedef MatrixType_ MatrixType;
144 typedef Preconditioner_ Preconditioner;
187template <
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
192 using Base::m_isInitialized;
193 using Base::m_iterations;
197 using Base::_solve_impl;
198 typedef MatrixType_ MatrixType;
199 typedef typename MatrixType::Scalar Scalar;
200 typedef typename MatrixType::RealScalar RealScalar;
201 typedef Preconditioner_ Preconditioner;
203 enum { UpLo = UpLo_ };
219 template <
typename MatrixDerived>
226 template <
typename Rhs,
typename Dest>
227 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
229 typedef typename Base::ActualMatrixType ActualMatrixType;
231 TransposeInput = (!MatrixWrapper::MatrixFree) && (UpLo == (
Lower |
Upper)) && (!MatrixType::IsRowMajor) &&
234 typedef std::conditional_t<TransposeInput, Transpose<const ActualMatrixType>, ActualMatrixType
const&>
236 EIGEN_STATIC_ASSERT(internal::check_implication(MatrixWrapper::MatrixFree, UpLo == (
Lower |
Upper)),
237 MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
238 typedef std::conditional_t<UpLo == (
Lower |
Upper), RowMajorWrapper,
239 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type>
243 m_error = Base::m_tolerance;
244 RowMajorWrapper row_mat(matrix());
245 internal::minres(SelfAdjointWrapper(row_mat), b, x, Base::m_preconditioner, m_iterations, m_error);
Index maxIterations() const
A minimal residual solver for sparse symmetric problems.
Definition MINRES.h:188
MINRES()
Definition MINRES.h:207
MINRES(const EigenBase< MatrixDerived > &A)
Definition MINRES.h:220
~MINRES()
Definition MINRES.h:223
Namespace containing all symbols from the Eigen library.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index