aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Eigenvalues
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2017-01-31 14:22:42 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2017-01-31 14:22:42 +0100
commit53026d29d41e81065b28631445e8eb5c4044c187 (patch)
tree1f38ee49f72cb17a89711fb0c4fe6160945bb986 /Eigen/src/Eigenvalues
parent63de19c0004933c7b2b1e418292b9f2ae6c138f4 (diff)
bug #478: fix regression in the eigen decomposition of zero matrices.
Diffstat (limited to 'Eigen/src/Eigenvalues')
-rw-r--r--Eigen/src/Eigenvalues/ComplexEigenSolver.h6
-rw-r--r--Eigen/src/Eigenvalues/RealSchur.h12
2 files changed, 16 insertions, 2 deletions
diff --git a/Eigen/src/Eigenvalues/ComplexEigenSolver.h b/Eigen/src/Eigenvalues/ComplexEigenSolver.h
index ec3b1633e..dc5fae06a 100644
--- a/Eigen/src/Eigenvalues/ComplexEigenSolver.h
+++ b/Eigen/src/Eigenvalues/ComplexEigenSolver.h
@@ -250,7 +250,7 @@ template<typename _MatrixType> class ComplexEigenSolver
EigenvectorType m_matX;
private:
- void doComputeEigenvectors(const RealScalar& matrixnorm);
+ void doComputeEigenvectors(RealScalar matrixnorm);
void sortEigenvalues(bool computeEigenvectors);
};
@@ -284,10 +284,12 @@ ComplexEigenSolver<MatrixType>::compute(const EigenBase<InputType>& matrix, bool
template<typename MatrixType>
-void ComplexEigenSolver<MatrixType>::doComputeEigenvectors(const RealScalar& matrixnorm)
+void ComplexEigenSolver<MatrixType>::doComputeEigenvectors(RealScalar matrixnorm)
{
const Index n = m_eivalues.size();
+ matrixnorm = numext::maxi(matrixnorm,(std::numeric_limits<RealScalar>::min)());
+
// Compute X such that T = X D X^(-1), where D is the diagonal of T.
// The matrix X is unit triangular.
m_matX = EigenvectorType::Zero(n, n);
diff --git a/Eigen/src/Eigenvalues/RealSchur.h b/Eigen/src/Eigenvalues/RealSchur.h
index d6a339f07..f5c86041d 100644
--- a/Eigen/src/Eigenvalues/RealSchur.h
+++ b/Eigen/src/Eigenvalues/RealSchur.h
@@ -248,12 +248,24 @@ template<typename MatrixType>
template<typename InputType>
RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU)
{
+ const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
+
eigen_assert(matrix.cols() == matrix.rows());
Index maxIters = m_maxIters;
if (maxIters == -1)
maxIters = m_maxIterationsPerRow * matrix.rows();
Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
+ if(scale<considerAsZero)
+ {
+ m_matT.setZero(matrix.rows(),matrix.cols());
+ if(computeU)
+ m_matU.setIdentity(matrix.rows(),matrix.cols());
+ m_info = Success;
+ m_isInitialized = true;
+ m_matUisUptodate = computeU;
+ return *this;
+ }
// Step 1. Reduce to Hessenberg form
m_hess.compute(matrix.derived()/scale);