Go to the documentation of this file.
55 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
69 if(rhs.norm() <= considerAsZero)
77 const Index maxIters = iters;
96 FMatrixType
H = FMatrixType::Zero(
m, restart + 1);
98 VectorType tau = VectorType::Zero(restart + 1);
101 std::vector < JacobiRotation < Scalar > >
G(restart);
109 r0.makeHouseholder(H0_tail, tau.coeffRef(0),
beta);
112 for (
Index k = 1; k <= restart; ++k)
116 v = VectorType::Unit(
m, k - 1);
120 for (
Index i = k - 1;
i >= 0; --
i) {
121 v.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
125 t.noalias() =
mat *
v;
126 v = precond.solve(
t);
131 v.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
134 if (
v.tail(
m - k).norm() != 0.0)
140 v.tail(
m - k).makeHouseholder(Hk_tail, tau.coeffRef(k),
beta);
143 v.tail(
m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
159 G[k - 1].makeGivens(
v(k - 1),
v(k));
162 v.applyOnTheLeft(k - 1, k,
G[k - 1].
adjoint());
163 w.applyOnTheLeft(k - 1, k,
G[k - 1].
adjoint());
167 H.col(k-1).head(k) =
v.head(k);
169 tol_error =
abs(
w(k)) / r0Norm;
170 bool stop = (k==
m || tol_error <
tol || iters == maxIters);
172 if (stop || k == restart)
176 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(
y);
184 x_new.tail(
m -
i).applyHouseholderOnTheLeft(
H.col(
i).tail(
m -
i - 1), tau.coeffRef(
i), workspace.data());
198 p0.noalias() = rhs -
mat*
x;
199 r0 = precond.solve(
p0);
207 r0.makeHouseholder(H0_tail, tau.coeffRef(0),
beta);
219 template<
typename _MatrixType,
220 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
225 template<
typename _MatrixType,
typename _Preconditioner>
268 template<
typename _MatrixType,
typename _Preconditioner>
303 template<
typename MatrixDerived>
318 template<
typename Rhs,
typename Dest>
335 #endif // EIGEN_GMRES_H
A GMRES solver for sparse square problems.
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate set rrange[*:*] noreverse nowriteback set trange[*:*] noreverse nowriteback set urange[*:*] noreverse nowriteback set vrange[*:*] noreverse nowriteback set xlabel matrix size set x2label set timefmt d m y n H
Namespace containing all symbols from the Eigen library.
void _solve_impl(const Rhs &b, Dest &x) const
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
double beta(double a, double b)
GMRES(const EigenBase< MatrixDerived > &A)
Index maxIterations() const
_Preconditioner Preconditioner
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, const Index &restart, typename Dest::RealScalar &tol_error)
_Preconditioner Preconditioner
void set_restart(const Index restart)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Base class for linear iterative solvers.
NumTraits< Scalar >::Real RealScalar
A matrix or vector expression mapping an existing expression.
IterativeSolverBase< GMRES > Base
JacobiRotation< float > G
const ActualMatrixType & matrix() const
Array< int, Dynamic, 1 > v
MatrixType::Scalar Scalar
Preconditioner m_preconditioner
void adjoint(const MatrixType &m)
Jet< T, N > sqrt(const Jet< T, N > &f)
MatrixType::RealScalar RealScalar
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:11:37