13 #ifndef EIGEN_MINRES_H
14 #define EIGEN_MINRES_H
17 #include "./InternalHeaderCheck.h"
32 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
34 void minres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
35 const Preconditioner& precond,
Index& iters,
36 typename Dest::RealScalar& tol_error)
39 typedef typename Dest::RealScalar RealScalar;
40 typedef typename Dest::Scalar Scalar;
44 const RealScalar rhsNorm2(rhs.squaredNorm());
54 const Index maxIters(iters);
55 const Index N(mat.cols());
56 const RealScalar threshold2(tol_error*tol_error*rhsNorm2);
60 VectorType v( VectorType::Zero(N) );
61 VectorType v_new(rhs-mat*x);
62 RealScalar residualNorm2(v_new.squaredNorm());
64 VectorType w_new(precond.solve(v_new));
66 RealScalar beta_new2(v_new.dot(w_new));
67 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
68 RealScalar beta_new(
sqrt(beta_new2));
69 const RealScalar beta_one(beta_new);
72 RealScalar c_old(1.0);
74 RealScalar s_old(0.0);
76 VectorType p_old(VectorType::Zero(N));
81 while ( iters < maxIters )
93 const RealScalar beta(beta_new);
99 v_new.noalias() = mat*w - beta*v_old;
100 const RealScalar alpha = v_new.dot(w);
102 w_new = precond.solve(v_new);
103 beta_new2 = v_new.dot(w_new);
104 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
105 beta_new =
sqrt(beta_new2);
108 const RealScalar r2 =s*alpha+c*c_old*beta;
109 const RealScalar r3 =s_old*beta;
110 const RealScalar r1_hat=c*alpha-c_old*s*beta;
111 const RealScalar r1 =
sqrt( std::pow(r1_hat,2) + std::pow(beta_new,2) );
120 p.noalias()=(w-r2*p_old-r3*p_oold) /r1;
121 x += beta_one*c*eta*p;
125 residualNorm2 *= s*s;
127 if ( residualNorm2 < threshold2)
138 tol_error = std::sqrt(residualNorm2 / rhsNorm2);
143 template<
typename MatrixType_,
int UpLo_=
Lower,
144 typename Preconditioner_ = IdentityPreconditioner>
149 template<
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
150 struct traits<MINRES<MatrixType_,UpLo_,Preconditioner_> >
152 typedef MatrixType_ MatrixType;
153 typedef Preconditioner_ Preconditioner;
196 template<
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
203 using Base::m_iterations;
205 using Base::m_isInitialized;
207 using Base::_solve_impl;
208 typedef MatrixType_ MatrixType;
209 typedef typename MatrixType::Scalar Scalar;
210 typedef typename MatrixType::RealScalar RealScalar;
211 typedef Preconditioner_ Preconditioner;
230 template<
typename MatrixDerived>
237 template<
typename Rhs,
typename Dest>
238 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
241 typedef typename Base::ActualMatrixType ActualMatrixType;
243 TransposeInput = (!MatrixWrapper::MatrixFree)
245 && (!MatrixType::IsRowMajor)
248 typedef std::conditional_t<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType
const&> RowMajorWrapper;
249 EIGEN_STATIC_ASSERT(internal::check_implication(MatrixWrapper::MatrixFree, UpLo==(
Lower|
Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
250 typedef std::conditional_t<UpLo==(
Lower|
Upper),
252 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
253 > SelfAdjointWrapper;
256 m_error = Base::m_tolerance;
257 RowMajorWrapper row_mat(matrix());
258 internal::minres(SelfAdjointWrapper(row_mat), b, x,
259 Base::m_preconditioner, m_iterations, m_error);
Index maxIterations() const
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:198
MINRES()
Definition: MINRES.h:218
MINRES(const EigenBase< MatrixDerived > &A)
Definition: MINRES.h:231
~MINRES()
Definition: MINRES.h:234
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)