57 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
59 Index &iters,
const Index &restart,
typename Dest::RealScalar & tol_error) {
64 typedef typename Dest::RealScalar RealScalar;
65 typedef typename Dest::Scalar Scalar;
71 if(rhs.norm() <= considerAsZero)
78 RealScalar tol = tol_error;
79 const Index maxIters = iters;
85 VectorType
p0 = rhs -
mat*
x;
86 VectorType r0 = precond.solve(
p0);
88 const RealScalar r0Norm = r0.norm();
98 FMatrixType
H = FMatrixType::Zero(
m, restart + 1);
99 VectorType
w = VectorType::Zero(restart + 1);
100 VectorType tau = VectorType::Zero(restart + 1);
103 std::vector < JacobiRotation < Scalar > >
G(restart);
106 VectorType t(
m),
v(
m), workspace(
m), x_new(
m);
111 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
114 for (
Index k = 1; k <= restart; ++k)
118 v = VectorType::Unit(
m, k - 1);
122 for (
Index i = k - 1;
i >= 0; --
i) {
123 v.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
127 t.noalias() =
mat *
v;
128 v = precond.solve(t);
133 v.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
136 if (
v.tail(
m - k).norm() != 0.0)
142 v.tail(
m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
145 v.tail(
m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
154 v.applyOnTheLeft(
i,
i + 1,
G[
i].adjoint());
158 if (k<
m &&
v(k) != (Scalar) 0)
161 G[k - 1].makeGivens(
v(k - 1),
v(k));
164 v.applyOnTheLeft(k - 1, k,
G[k - 1].adjoint());
165 w.applyOnTheLeft(k - 1, k,
G[k - 1].adjoint());
169 H.col(k-1).head(k) =
v.head(k);
171 tol_error =
abs(
w(k)) / r0Norm;
172 bool stop = (k==
m || tol_error < tol || iters == maxIters);
174 if (stop || k == restart)
178 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(
y);
186 x_new.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
200 p0.noalias() = rhs -
mat*
x;
201 r0 = precond.solve(
p0);
209 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
221 template<
typename MatrixType_,
222 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
227 template<
typename MatrixType_,
typename Preconditioner_>
228 struct traits<GMRES<MatrixType_,Preconditioner_> >
231 typedef Preconditioner_ Preconditioner;
270 template<
typename MatrixType_,
typename Preconditioner_>
286 typedef typename MatrixType::Scalar
Scalar;
305 template<
typename MatrixDerived>
320 template<
typename Rhs,
typename Dest>
Array< int, Dynamic, 1 > v
SparseMatrix< double > A(n, n)
JacobiRotation< float > G
Matrix< float, 1, Dynamic > MatrixType
A GMRES solver for sparse square problems.
MatrixType::RealScalar RealScalar
GMRES(const EigenBase< MatrixDerived > &A)
const ActualMatrixType & matrix() const
Preconditioner_ Preconditioner
IterativeSolverBase< GMRES > Base
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
MatrixType::Scalar Scalar
void set_restart(const Index restart)
Index maxIterations() const
void _solve_impl(const Rhs &b, Dest &x) const
Preconditioner m_preconditioner
const ActualMatrixType & matrix() const
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, const Index &restart, 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
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)