mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-05-01 16:24:28 +08:00
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.
This commit is contained in:
parent
d655900953
commit
9f0c5c3669
@ -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();
|
||||
|
Loading…
x
Reference in New Issue
Block a user