13 #ifndef EIGEN_MINRES_H
14 #define EIGEN_MINRES_H
32 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
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);
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;
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_> >
153 typedef Preconditioner_ Preconditioner;
196 template<
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
209 typedef typename MatrixType::Scalar
Scalar;
230 template<
typename MatrixDerived>
237 template<
typename Rhs,
typename Dest>
243 TransposeInput = (!MatrixWrapper::MatrixFree)
245 && (!MatrixType::IsRowMajor)
248 typedef std::conditional_t<TransposeInput,Transpose<const ActualMatrixType>,
ActualMatrixType const&> RowMajorWrapper;
252 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
253 > SelfAdjointWrapper;
257 RowMajorWrapper row_mat(
matrix());
Array< int, Dynamic, 1 > v
SparseMatrix< double > A(n, n)
#define EIGEN_DONT_INLINE
#define EIGEN_STATIC_ASSERT(X, MSG)
Matrix< float, 1, Dynamic > MatrixType
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
Index maxIterations() const
MatrixWrapper::ActualMatrixType ActualMatrixType
void _solve_impl(const Rhs &b, Dest &x) const
Preconditioner m_preconditioner
const ActualMatrixType & matrix() const
A minimal residual solver for sparse symmetric problems.
const ActualMatrixType & matrix() const
Preconditioner_ Preconditioner
MatrixType::RealScalar RealScalar
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
MINRES(const EigenBase< MatrixDerived > &A)
MatrixType::Scalar Scalar
IterativeSolverBase< MINRES > Base
constexpr bool check_implication(bool a, bool b)
EIGEN_DONT_INLINE void minres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Eigen::internal::remove_all_t< DerType >, typename Eigen::internal::traits< Eigen::internal::remove_all_t< DerType >>::Scalar, product) > sqrt(const Eigen::AutoDiffScalar< DerType > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(internal::remove_all_t< DerType >, typename internal::traits< internal::remove_all_t< DerType >>::Scalar, product) > pow(const Eigen::AutoDiffScalar< DerType > &x, const typename internal::traits< internal::remove_all_t< DerType >>::Scalar &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)