aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseQR
diff options
context:
space:
mode:
authorGravatar Jeff Trull <edaskel@att.net>2018-02-15 15:00:31 -0800
committerGravatar Jeff Trull <edaskel@att.net>2018-02-15 15:00:31 -0800
commit9f0c5c3669a985a0dfa44db2ab05306edf1130c1 (patch)
tree93812560d0404703609acf173504e6bca8c22ea0 /Eigen/src/SparseQR
parentd6559009530b87f80fc5aaa864c012ae5109e848 (diff)
Make sparse QR result sizes consistent with dense QR, with the following rules:
1) Q is always square 2) Q*R*P' is valid and recovers the original matrix This implies that the size of Q is the number of rows in the original matrix, square, and that the size of R is the size of the original matrix.
Diffstat (limited to 'Eigen/src/SparseQR')
-rw-r--r--Eigen/src/SparseQR/SparseQR.h13
1 files changed, 8 insertions, 5 deletions
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index 278d775aa..7409fcae9 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -328,7 +328,7 @@ void SparseQR<MatrixType,OrderingType>::analyzePattern(const MatrixType& mat)
internal::coletree(matCpy, m_etree, m_firstRowElt, m_outputPerm_c.indices().data());
m_isEtreeOk = true;
- m_R.resize(diagSize, n);
+ m_R.resize(m, n);
m_Q.resize(m, diagSize);
// Allocate space for nonzero elements : rough estimation
@@ -605,7 +605,7 @@ struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived
// Get the references
SparseQR_QProduct(const SparseQRType& qr, const Derived& other, bool transpose) :
m_qr(qr),m_other(other),m_transpose(transpose) {}
- inline Index rows() const { return m_transpose ? m_qr.rows() : m_qr.cols(); }
+ inline Index rows() const { return m_qr.matrixQ().rows(); }
inline Index cols() const { return m_other.cols(); }
// Assign to a vector
@@ -633,7 +633,10 @@ struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived
}
else
{
- eigen_assert(m_qr.m_Q.rows() == m_other.rows() && "Non conforming object sizes");
+ eigen_assert(m_qr.matrixQ().cols() == m_other.rows() && "Non conforming object sizes");
+
+ res.conservativeResize(rows(), cols());
+
// Compute res = Q * other column by column
for(Index j = 0; j < res.cols(); j++)
{
@@ -675,7 +678,7 @@ struct SparseQRMatrixQReturnType : public EigenBase<SparseQRMatrixQReturnType<Sp
return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
}
inline Index rows() const { return m_qr.rows(); }
- inline Index cols() const { return (std::min)(m_qr.rows(),m_qr.cols()); }
+ inline Index cols() const { return m_qr.rows(); }
// To use for operations with the transpose of Q FIXME this is the same as adjoint at the moment
SparseQRMatrixQTransposeReturnType<SparseQRType> transpose() const
{
@@ -715,7 +718,7 @@ struct Assignment<DstXprType, SparseQRMatrixQReturnType<SparseQRType>, internal:
typedef typename DstXprType::StorageIndex StorageIndex;
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,Scalar> &/*func*/)
{
- typename DstXprType::PlainObject idMat(src.m_qr.rows(), src.m_qr.rows());
+ typename DstXprType::PlainObject idMat(src.rows(), src.cols());
idMat.setIdentity();
// Sort the sparse householder reflectors if needed
const_cast<SparseQRType *>(&src.m_qr)->_sort_matrix_Q();