aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseQR
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2013-09-12 22:16:35 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2013-09-12 22:16:35 +0200
commit1b4623e7132c0959925386801e9e342ba5565f19 (patch)
treeda5c44fba65ad414ff667b8604b66dbab6c2a5d1 /Eigen/src/SparseQR
parent4612a1cd870b4a72b3849608fcdce9a18dc80a65 (diff)
Fix elimination tree and SparseQR with rows<cols
Diffstat (limited to 'Eigen/src/SparseQR')
-rw-r--r--Eigen/src/SparseQR/SparseQR.h9
1 files changed, 5 insertions, 4 deletions
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index c42dfc556..50c4c1145 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -161,8 +161,9 @@ class SparseQR
b = y;
// Solve with the triangular matrix R
+ y.resize((std::max)(cols(),Index(y.rows())),y.cols());
y.topRows(rank) = this->matrixR().topLeftCorner(rank, rank).template triangularView<Upper>().solve(b.topRows(rank));
- y.bottomRows(y.size()-rank).setZero();
+ y.bottomRows(y.rows()-rank).setZero();
// Apply the column permutation
if (m_perm_c.size()) dest.topRows(cols()) = colsPermutation() * y.topRows(cols());
@@ -246,7 +247,7 @@ class SparseQR
Index m_nonzeropivots; // Number of non zero pivots found
IndexVector m_etree; // Column elimination tree
IndexVector m_firstRowElt; // First element in each row
- bool m_isQSorted; // whether Q is sorted or not
+ bool m_isQSorted; // whether Q is sorted or not
template <typename, typename > friend struct SparseQR_QProduct;
template <typename > friend struct SparseQRMatrixQReturnType;
@@ -338,7 +339,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
Index nonzeroCol = 0; // Record the number of valid pivots
// 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 < (std::min)(n,m); ++col)
{
mark.setConstant(-1);
m_R.startVec(col);
@@ -346,7 +347,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
mark(nonzeroCol) = col;
Qidx(0) = nonzeroCol;
nzcolR = 0; nzcolQ = 1;
- found_diag = false;
+ found_diag = col>=m;
tval.setZero();
// Symbolic factorization: find the nonzero locations of the column k of the factors R and Q, i.e.,