aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseQR
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-10-19 22:42:20 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-10-19 22:42:20 +0200
commit8838b0a1fff2cb01fab3a55312a16ae20daead13 (patch)
tree8eff174e86bf31105a79a33957451e3b06f10b6c /Eigen/src/SparseQR
parentb50e5bc816905e27423d64a9685c9015016d6da5 (diff)
Fix SparseQR::rank for a completely empty matrix.
Diffstat (limited to 'Eigen/src/SparseQR')
-rw-r--r--Eigen/src/SparseQR/SparseQR.h6
1 files changed, 4 insertions, 2 deletions
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index f6a9af672..879ac553d 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -378,6 +378,8 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
{
RealScalar max2Norm = 0.0;
for (int j = 0; j < n; j++) max2Norm = (max)(max2Norm, m_pmat.col(j).norm());
+ if(max2Norm==RealScalar(0))
+ max2Norm = RealScalar(1);
pivotThreshold = 20 * (m + n) * max2Norm * NumTraits<RealScalar>::epsilon();
}
@@ -386,7 +388,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
Index nonzeroCol = 0; // Record the number of valid pivots
m_Q.startVec(0);
-
+
// Left looking rank-revealing QR factorization: compute a column of R and Q at a time
for (Index col = 0; col < n; ++col)
{
@@ -554,7 +556,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
m_R.finalize();
m_R.makeCompressed();
m_isQSorted = false;
-
+
m_nonzeropivots = nonzeroCol;
if(nonzeroCol<n)