10 #ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H 11 #define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H 16 template <
typename _MatrixType>
47 template <
typename _MatrixType>
52 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
53 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
54 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
55 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
94 : m_cpqr(rows, cols), m_zCoeffs((
std::
min)(rows, cols)), m_temp(cols) {}
112 template <
typename InputType>
114 : m_cpqr(matrix.
rows(), matrix.
cols()),
116 m_temp(matrix.
cols())
127 template<
typename InputType>
129 : m_cpqr(matrix.derived()),
131 m_temp(matrix.
cols())
146 template <
typename Rhs>
150 "CompleteOrthogonalDecomposition is not initialized.");
160 MatrixType
Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
161 applyZAdjointOnTheLeftInPlace(Z);
168 const MatrixType&
matrixQTZ()
const {
return m_cpqr.matrixQR(); }
181 const MatrixType&
matrixT()
const {
return m_cpqr.matrixQR(); }
183 template <
typename InputType>
193 return m_cpqr.colsPermutation();
289 inline const HCoeffsType&
hCoeffs()
const {
return m_cpqr.hCoeffs(); }
296 const HCoeffsType&
zCoeffs()
const {
return m_zCoeffs; }
339 RealScalar
threshold()
const {
return m_cpqr.threshold(); }
353 inline RealScalar
maxPivot()
const {
return m_cpqr.maxPivot(); }
364 eigen_assert(m_cpqr.m_isInitialized &&
"Decomposition is not initialized.");
368 #ifndef EIGEN_PARSED_BY_DOXYGEN 369 template <
typename RhsType,
typename DstType>
370 EIGEN_DEVICE_FUNC
void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
378 void computeInPlace();
382 template <
typename Rhs>
383 void applyZAdjointOnTheLeftInPlace(
Rhs& rhs)
const;
390 template <
typename MatrixType>
393 return m_cpqr.absDeterminant();
396 template <
typename MatrixType>
399 return m_cpqr.logAbsDeterminant();
409 template <
typename MatrixType>
412 check_template_parameters();
417 const Index rank = m_cpqr.rank();
420 m_zCoeffs.resize((
std::min)(rows, cols));
435 for (
Index k = rank - 1; k >= 0; --k) {
440 m_cpqr.m_qr.col(k).head(k + 1).swap(
441 m_cpqr.m_qr.col(rank - 1).head(k + 1));
448 .tail(cols - rank + 1)
449 .makeHouseholderInPlace(m_zCoeffs(k), beta);
450 m_cpqr.m_qr(k, rank - 1) = beta;
453 m_cpqr.m_qr.topRightCorner(k, cols - rank + 1)
454 .applyHouseholderOnTheRight(
455 m_cpqr.m_qr.row(k).tail(cols - rank).transpose(), m_zCoeffs(k),
460 m_cpqr.m_qr.col(k).head(k + 1).swap(
461 m_cpqr.m_qr.col(rank - 1).head(k + 1));
467 template <
typename MatrixType>
468 template <
typename Rhs>
472 const Index nrhs = rhs.cols();
473 const Index rank = this->rank();
475 for (
Index k = 0; k < rank; ++k) {
477 rhs.row(k).swap(rhs.row(rank - 1));
479 rhs.middleRows(rank - 1, cols - rank + 1)
480 .applyHouseholderOnTheLeft(
484 rhs.row(k).swap(rhs.row(rank - 1));
489 #ifndef EIGEN_PARSED_BY_DOXYGEN 490 template <
typename _MatrixType>
491 template <
typename RhsType,
typename DstType>
493 const RhsType& rhs, DstType& dst)
const {
496 const Index rank = this->rank();
505 typename RhsType::PlainObject
c(rhs);
510 dst.topRows(rank) = matrixT()
511 .topLeftCorner(rank, rank)
512 .template triangularView<Upper>()
513 .solve(c.topRows(rank));
519 dst.bottomRows(cols - rank).setZero();
520 applyZAdjointOnTheLeftInPlace(dst);
524 dst = colsPermutation() * dst;
530 template<
typename DstXprType,
typename MatrixType>
544 template <
typename MatrixType>
547 return m_cpqr.householderQ();
554 template <
typename Derived>
562 #endif // EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H MatrixType matrixZ() const
const HCoeffsType & zCoeffs() const
const Solve< CompleteOrthogonalDecomposition, Rhs > solve(const MatrixBase< Rhs > &b) const
HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
Inverse< CodType > SrcXprType
CompleteOrthogonalDecomposition()
Default Constructor.
EIGEN_DEVICE_FUNC Index rows() const
MatrixType::StorageIndex StorageIndex
void adjoint(const MatrixType &m)
EIGEN_DEVICE_FUNC SegmentReturnType tail(Index n)
This is the const version of tail(Index).
MatrixType::RealScalar absDeterminant() const
MatrixType::PlainObject PlainObject
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was succesful.
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Namespace containing all symbols from the Eigen library.
MatrixType::RealScalar logAbsDeterminant() const
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
const MatrixType & matrixQTZ() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Index dimensionOfKernel() const
bool isSurjective() const
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
const MatrixType & matrixT() const
HouseholderSequenceType householderQ(void) const
Complete orthogonal decomposition (COD) of a matrix.
Eigen::Index Index
The interface type of indices.
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
MatrixType::Scalar Scalar
MatrixType::RealScalar RealScalar
internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
internal::plain_diag_type< MatrixType >::type HCoeffsType
Sequence of Householder reflections acting on subspaces with decreasing size.
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
Expression of the inverse of another expression.
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
RealScalar maxPivot() const
HouseholderSequenceType matrixQ(void) const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
ColPivHouseholderQR< MatrixType > m_cpqr
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
const CompleteOrthogonalDecomposition< PlainObject > completeOrthogonalDecomposition() const
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationType
NumTraits< Scalar >::Real RealScalar
CompleteOrthogonalDecomposition & compute(const EigenBase< InputType > &matrix)
const PermutationType & colsPermutation() const
PermutationType::Index PermIndexType
internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
bool isInvertible() const
internal::nested_eval< T, 1 >::type eval(const T &xpr)
internal::conditional< NumTraits< Scalar >::IsComplex, const CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const Derived >, const Derived & >::type ConjugateReturnType
static void check_template_parameters()
Pseudo expression representing a solving operation.
CompleteOrthogonalDecomposition & setThreshold(Default_t)
Index nonzeroPivots() const
const HCoeffsType & hCoeffs() const
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
internal::plain_row_type< MatrixType >::type RowVectorType
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
The matrix class, also used for vectors and row-vectors.
EIGEN_DEVICE_FUNC const XprTypeNestedCleaned & nestedExpression() const
RealScalar threshold() const
EIGEN_DEVICE_FUNC Derived & derived()
Base class for all dense matrices, vectors, and expressions.
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< typename DstXprType::Scalar, typename CodType::Scalar > &)
CompleteOrthogonalDecomposition< MatrixType > CodType
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const