mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-19 16:19:37 +08:00
make ColPivotingQR use HouseholderSequence
This commit is contained in:
parent
49dd5d7847
commit
24950bdfcb
@ -77,11 +77,12 @@ template<typename VectorType, int Size, int PacketAccess> class VectorBlock
|
||||
typedef Block<VectorType,
|
||||
ei_traits<VectorType>::RowsAtCompileTime==1 ? 1 : Size,
|
||||
ei_traits<VectorType>::ColsAtCompileTime==1 ? 1 : Size,
|
||||
PacketAccess> Base;
|
||||
PacketAccess> _Base;
|
||||
enum {
|
||||
IsColVector = ei_traits<VectorType>::ColsAtCompileTime==1
|
||||
};
|
||||
public:
|
||||
_EIGEN_GENERIC_PUBLIC_INTERFACE(VectorBlock, _Base)
|
||||
|
||||
using Base::operator=;
|
||||
using Base::operator+=;
|
||||
|
@ -136,6 +136,22 @@ template<typename VectorsType, typename CoeffsType> class HouseholderSequence
|
||||
}
|
||||
}
|
||||
|
||||
template<typename OtherDerived>
|
||||
typename OtherDerived::PlainMatrixType operator*(const MatrixBase<OtherDerived>& other) const
|
||||
{
|
||||
typename OtherDerived::PlainMatrixType res(other);
|
||||
applyThisOnTheLeft(res);
|
||||
return res;
|
||||
}
|
||||
|
||||
template<typename OtherDerived> friend
|
||||
typename OtherDerived::PlainMatrixType operator*(const MatrixBase<OtherDerived>& other, const HouseholderSequence& h)
|
||||
{
|
||||
typename OtherDerived::PlainMatrixType res(other);
|
||||
h.applyThisOnTheRight(res);
|
||||
return res;
|
||||
}
|
||||
|
||||
protected:
|
||||
|
||||
typename VectorsType::Nested m_vectors;
|
||||
@ -143,4 +159,10 @@ template<typename VectorsType, typename CoeffsType> class HouseholderSequence
|
||||
bool m_trans;
|
||||
};
|
||||
|
||||
template<typename VectorsType, typename CoeffsType>
|
||||
HouseholderSequence<VectorsType,CoeffsType> makeHouseholderSequence(const VectorsType& v, const CoeffsType& h, bool trans=false)
|
||||
{
|
||||
return HouseholderSequence<VectorsType,CoeffsType>(v, h, trans);
|
||||
}
|
||||
|
||||
#endif // EIGEN_HOUSEHOLDER_SEQUENCE_H
|
||||
|
@ -62,6 +62,7 @@ template<typename MatrixType> class ColPivotingHouseholderQR
|
||||
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
|
||||
typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
|
||||
typedef Matrix<RealScalar, 1, ColsAtCompileTime> RealRowVectorType;
|
||||
typedef typename HouseholderSequence<MatrixQType,HCoeffsType>::ConjugateReturnType HouseholderSequenceType;
|
||||
|
||||
/**
|
||||
* \brief Default Constructor.
|
||||
@ -99,7 +100,7 @@ template<typename MatrixType> class ColPivotingHouseholderQR
|
||||
template<typename OtherDerived, typename ResultType>
|
||||
bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
|
||||
|
||||
MatrixQType matrixQ(void) const;
|
||||
HouseholderSequenceType matrixQ(void) const;
|
||||
|
||||
/** \returns a reference to the matrix where the Householder QR decomposition is stored
|
||||
*/
|
||||
@ -355,13 +356,8 @@ bool ColPivotingHouseholderQR<MatrixType>::solve(
|
||||
|
||||
typename OtherDerived::PlainMatrixType c(b);
|
||||
|
||||
Matrix<Scalar,1,MatrixType::ColsAtCompileTime> temp(cols);
|
||||
for (int k = 0; k < m_rank; ++k)
|
||||
{
|
||||
int remainingSize = rows-k;
|
||||
c.corner(BottomRight, remainingSize, cols)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(remainingSize-1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
|
||||
}
|
||||
// Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
|
||||
c.applyOnTheLeft(makeHouseholderSequence(m_qr.corner(TopLeft,rows,m_rank), m_hCoeffs.start(m_rank)).transpose());
|
||||
|
||||
if(!isSurjective())
|
||||
{
|
||||
@ -381,25 +377,12 @@ bool ColPivotingHouseholderQR<MatrixType>::solve(
|
||||
return true;
|
||||
}
|
||||
|
||||
/** \returns the matrix Q */
|
||||
/** \returns the matrix Q as a sequence of householder transformations */
|
||||
template<typename MatrixType>
|
||||
typename ColPivotingHouseholderQR<MatrixType>::MatrixQType ColPivotingHouseholderQR<MatrixType>::matrixQ() const
|
||||
typename ColPivotingHouseholderQR<MatrixType>::HouseholderSequenceType ColPivotingHouseholderQR<MatrixType>::matrixQ() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
|
||||
// compute the product H'_0 H'_1 ... H'_n-1,
|
||||
// where H_k is the k-th Householder transformation I - h_k v_k v_k'
|
||||
// and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
|
||||
int rows = m_qr.rows();
|
||||
int cols = m_qr.cols();
|
||||
int size = std::min(rows,cols);
|
||||
MatrixQType res = MatrixQType::Identity(rows, rows);
|
||||
Matrix<Scalar,1,MatrixType::RowsAtCompileTime> temp(rows);
|
||||
for (int k = size-1; k >= 0; k--)
|
||||
{
|
||||
res.block(k, k, rows-k, rows-k)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(rows-k-1), ei_conj(m_hCoeffs.coeff(k)), &temp.coeffRef(k));
|
||||
}
|
||||
return res;
|
||||
return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
|
||||
}
|
||||
|
||||
#endif // EIGEN_HIDE_HEAVY_CODE
|
||||
|
Loading…
x
Reference in New Issue
Block a user