mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-10 02:39:03 +08:00
Correct GMRES:
* Fix error in calculation of residual at restart. * Use relative residual as stopping criterion. * Improve documentation.
This commit is contained in:
parent
e51da9c3a8
commit
953ec08089
@ -27,7 +27,7 @@ namespace internal {
|
||||
* \param iters on input: maximum number of iterations to perform
|
||||
* on output: number of iterations performed
|
||||
* \param restart number of iterations for a restart
|
||||
* \param tol_error on input: residual tolerance
|
||||
* \param tol_error on input: relative residual tolerance
|
||||
* on output: residuum achieved
|
||||
*
|
||||
* \sa IterativeMethods::bicgstab()
|
||||
@ -70,18 +70,24 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition
|
||||
|
||||
const int m = mat.rows();
|
||||
|
||||
VectorType p0 = rhs - mat*x;
|
||||
// residual and preconditioned residual
|
||||
const VectorType p0 = rhs - mat*x;
|
||||
VectorType r0 = precond.solve(p0);
|
||||
|
||||
const RealScalar r0Norm = r0.norm();
|
||||
|
||||
// is initial guess already good enough?
|
||||
if(abs(r0.norm()) < tol) {
|
||||
if(r0Norm == 0) {
|
||||
tol_error=0;
|
||||
return true;
|
||||
}
|
||||
|
||||
// storage for Hessenberg matrix and Householder data
|
||||
FMatrixType H = FMatrixType::Zero(m, restart + 1);
|
||||
VectorType w = VectorType::Zero(restart + 1);
|
||||
|
||||
FMatrixType H = FMatrixType::Zero(m, restart + 1); // Hessenberg matrix
|
||||
VectorType tau = VectorType::Zero(restart + 1);
|
||||
|
||||
// storage for Jacobi rotations
|
||||
std::vector < JacobiRotation < Scalar > > G(restart);
|
||||
|
||||
// generate first Householder vector
|
||||
@ -112,7 +118,6 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition
|
||||
}
|
||||
|
||||
if (v.tail(m - k).norm() != 0.0) {
|
||||
|
||||
if (k <= restart) {
|
||||
|
||||
// generate new Householder vector
|
||||
@ -135,19 +140,19 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition
|
||||
}
|
||||
|
||||
if (k<m && v(k) != (Scalar) 0) {
|
||||
|
||||
// determine next Givens rotation
|
||||
G[k - 1].makeGivens(v(k - 1), v(k));
|
||||
|
||||
// apply Givens rotation to v and w
|
||||
v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
|
||||
w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
|
||||
|
||||
}
|
||||
|
||||
// insert coefficients into upper matrix triangle
|
||||
H.col(k - 1).head(k) = v.head(k);
|
||||
|
||||
bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
|
||||
bool stop=(k==m || abs(w(k)) < tol * r0Norm || iters == maxIters);
|
||||
|
||||
if (stop || k == restart) {
|
||||
|
||||
@ -174,13 +179,13 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition
|
||||
} else {
|
||||
k=0;
|
||||
|
||||
// reset data for a restart r0 = rhs - mat * x;
|
||||
VectorType p0=mat*x;
|
||||
VectorType p1=precond.solve(p0);
|
||||
r0 = rhs - p1;
|
||||
// r0_sqnorm = r0.squaredNorm();
|
||||
w = VectorType::Zero(restart + 1);
|
||||
// reset data for restart
|
||||
const VectorType p0 = rhs - mat*x;
|
||||
r0 = precond.solve(p0);
|
||||
|
||||
// clear Hessenberg matrix and Householder data
|
||||
H = FMatrixType::Zero(m, restart + 1);
|
||||
w = VectorType::Zero(restart + 1);
|
||||
tau = VectorType::Zero(restart + 1);
|
||||
|
||||
// generate first Householder vector
|
||||
@ -194,7 +199,6 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition
|
||||
}
|
||||
|
||||
|
||||
|
||||
}
|
||||
|
||||
return false;
|
||||
|
Loading…
x
Reference in New Issue
Block a user