diff options
author | Kolja Brix <brix@igpm.rwth-aachen.de> | 2014-07-10 08:20:55 +0200 |
---|---|---|
committer | Kolja Brix <brix@igpm.rwth-aachen.de> | 2014-07-10 08:20:55 +0200 |
commit | e955725ff1434396cc41beda6f5989bef0940704 (patch) | |
tree | 4422efaf1f6360f5b83f3ca4b6d883f410fc958c /unsupported/Eigen/src | |
parent | 23bb592a2d0de937c68fc1c819642bbcfa2c1e4f (diff) |
Fix GMRES: Initialize essential Householder vector with correct dimension. Add check if initial guess is already a sufficient approximation.
Diffstat (limited to 'unsupported/Eigen/src')
-rw-r--r-- | unsupported/Eigen/src/IterativeSolvers/GMRES.h | 12 |
1 files changed, 8 insertions, 4 deletions
diff --git a/unsupported/Eigen/src/IterativeSolvers/GMRES.h b/unsupported/Eigen/src/IterativeSolvers/GMRES.h index 073367506..c8c84069e 100644 --- a/unsupported/Eigen/src/IterativeSolvers/GMRES.h +++ b/unsupported/Eigen/src/IterativeSolvers/GMRES.h @@ -2,7 +2,7 @@ // for linear algebra. // // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr> -// Copyright (C) 2012 Kolja Brix <brix@igpm.rwth-aaachen.de> +// Copyright (C) 2012, 2014 Kolja Brix <brix@igpm.rwth-aaachen.de> // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed @@ -72,16 +72,20 @@ bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Precondition VectorType p0 = rhs - mat*x; VectorType r0 = precond.solve(p0); -// RealScalar r0_sqnorm = r0.squaredNorm(); + + // is initial guess already good enough? + if(abs(r0.norm()) < tol) { + return true; + } VectorType w = VectorType::Zero(restart + 1); - FMatrixType H = FMatrixType::Zero(m, restart + 1); + FMatrixType H = FMatrixType::Zero(m, restart + 1); // Hessenberg matrix VectorType tau = VectorType::Zero(restart + 1); std::vector < JacobiRotation < Scalar > > G(restart); // generate first Householder vector - VectorType e; + VectorType e(m-1); RealScalar beta; r0.makeHouseholder(e, tau.coeffRef(0), beta); w(0)=(Scalar) beta; |