mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-09-12 17:33:15 +08:00
Fix SparseQR::rank for a completely empty matrix.
This commit is contained in:
parent
b50e5bc816
commit
8838b0a1ff
@ -378,6 +378,8 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
|
|||||||
{
|
{
|
||||||
RealScalar max2Norm = 0.0;
|
RealScalar max2Norm = 0.0;
|
||||||
for (int j = 0; j < n; j++) max2Norm = (max)(max2Norm, m_pmat.col(j).norm());
|
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();
|
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
|
Index nonzeroCol = 0; // Record the number of valid pivots
|
||||||
m_Q.startVec(0);
|
m_Q.startVec(0);
|
||||||
|
|
||||||
// Left looking rank-revealing QR factorization: compute a column of R and Q at a time
|
// Left looking rank-revealing QR factorization: compute a column of R and Q at a time
|
||||||
for (Index col = 0; col < n; ++col)
|
for (Index col = 0; col < n; ++col)
|
||||||
{
|
{
|
||||||
@ -554,7 +556,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
|
|||||||
m_R.finalize();
|
m_R.finalize();
|
||||||
m_R.makeCompressed();
|
m_R.makeCompressed();
|
||||||
m_isQSorted = false;
|
m_isQSorted = false;
|
||||||
|
|
||||||
m_nonzeropivots = nonzeroCol;
|
m_nonzeropivots = nonzeroCol;
|
||||||
|
|
||||||
if(nonzeroCol<n)
|
if(nonzeroCol<n)
|
||||||
|
Loading…
x
Reference in New Issue
Block a user