10#ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
11#define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
14#include "./InternalHeaderCheck.h"
19template <
typename MatrixType_,
typename PermutationIndex_>
20struct traits<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_>> : traits<MatrixType_> {
21 typedef MatrixXpr XprKind;
22 typedef SolverStorage StorageKind;
23 typedef PermutationIndex_ PermutationIndex;
52template <
typename MatrixType_,
typename PermutationIndex_>
54 :
public SolverBase<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_>> {
56 typedef MatrixType_ MatrixType;
59 template <
typename Derived>
60 friend struct internal::solve_assertion;
61 typedef PermutationIndex_ PermutationIndex;
64 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
65 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
67 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
69 typedef typename internal::plain_row_type<MatrixType, Index>::type IntRowVectorType;
70 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
71 typedef typename internal::plain_row_type<MatrixType, RealScalar>::type RealRowVectorType;
73 HouseholderSequenceType;
74 typedef typename MatrixType::PlainObject PlainObject;
93 : m_cpqr(rows, cols), m_zCoeffs((std::min)(rows, cols)), m_temp(cols) {}
111 template <
typename InputType>
113 : m_cpqr(matrix.rows(), matrix.cols()),
114 m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
115 m_temp(matrix.cols()) {
126 template <
typename InputType>
128 : m_cpqr(matrix.
derived()), m_zCoeffs((std::min)(matrix.rows(), matrix.cols())), m_temp(matrix.cols()) {
132#ifdef EIGEN_PARSED_BY_DOXYGEN
142 template <
typename Rhs>
147 HouseholderSequenceType matrixQ(
void)
const {
return m_cpqr.householderQ(); }
152 MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
160 const MatrixType&
matrixQTZ()
const {
return m_cpqr.matrixQR(); }
173 const MatrixType&
matrixT()
const {
return m_cpqr.matrixQR(); }
175 template <
typename InputType>
178 m_cpqr.compute(matrix);
297 eigen_assert(m_cpqr.m_isInitialized &&
"CompleteOrthogonalDecomposition is not initialized.");
301 inline Index rows()
const {
return m_cpqr.rows(); }
302 inline Index cols()
const {
return m_cpqr.cols(); }
309 inline const HCoeffsType&
hCoeffs()
const {
return m_cpqr.hCoeffs(); }
316 const HCoeffsType&
zCoeffs()
const {
return m_zCoeffs; }
359 RealScalar
threshold()
const {
return m_cpqr.threshold(); }
373 inline RealScalar
maxPivot()
const {
return m_cpqr.maxPivot(); }
384 eigen_assert(m_cpqr.m_isInitialized &&
"Decomposition is not initialized.");
388#ifndef EIGEN_PARSED_BY_DOXYGEN
389 template <
typename RhsType,
typename DstType>
390 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
392 template <
bool Conjugate,
typename RhsType,
typename DstType>
393 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
397 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar)
399 template <
bool Transpose_,
typename Rhs>
400 void _check_solve_assertion(
const Rhs& b)
const {
401 EIGEN_ONLY_USED_FOR_DEBUG(b);
402 eigen_assert(m_cpqr.m_isInitialized &&
"CompleteOrthogonalDecomposition is not initialized.");
403 eigen_assert((Transpose_ ?
derived().cols() :
derived().rows()) == b.rows() &&
404 "CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
413 template <
bool Conjugate,
typename Rhs>
418 template <
typename Rhs>
422 HCoeffsType m_zCoeffs;
423 RowVectorType m_temp;
426template <
typename MatrixType,
typename PermutationIndex>
431template <
typename MatrixType,
typename PermutationIndex>
433 return m_cpqr.absDeterminant();
436template <
typename MatrixType,
typename PermutationIndex>
439 return m_cpqr.logAbsDeterminant();
442template <
typename MatrixType,
typename PermutationIndex>
444 return m_cpqr.signDeterminant();
454template <
typename MatrixType,
typename PermutationIndex>
456 eigen_assert(m_cpqr.cols() <= NumTraits<PermutationIndex>::highest());
459 const Index cols = m_cpqr.cols();
460 const Index rows = m_cpqr.rows();
461 m_zCoeffs.resize((std::min)(rows, cols));
481 m_cpqr.m_qr.col(k).head(k + 1).swap(m_cpqr.m_qr.col(
rank - 1).head(k + 1));
487 m_cpqr.m_qr.row(k).tail(cols -
rank + 1).makeHouseholderInPlace(m_zCoeffs(k), beta);
488 m_cpqr.m_qr(k,
rank - 1) = beta;
491 m_cpqr.m_qr.topRightCorner(k, cols -
rank + 1)
492 .applyHouseholderOnTheRight(m_cpqr.m_qr.row(k).tail(cols -
rank).adjoint(), m_zCoeffs(k), &m_temp(0));
496 m_cpqr.m_qr.col(k).head(k + 1).swap(m_cpqr.m_qr.col(
rank - 1).head(k + 1));
502template <
typename MatrixType,
typename PermutationIndex>
503template <
bool Conjugate,
typename Rhs>
505 const Index cols = this->cols();
506 const Index nrhs = rhs.cols();
511 rhs.row(k).swap(rhs.row(
rank - 1));
513 rhs.middleRows(
rank - 1, cols -
rank + 1)
514 .applyHouseholderOnTheLeft(
matrixQTZ().row(k).tail(cols -
rank).
transpose().
template conjugateIf<!Conjugate>(),
515 zCoeffs().
template conjugateIf<Conjugate>()(k), &temp(0));
517 rhs.row(k).swap(rhs.row(
rank - 1));
522template <
typename MatrixType,
typename PermutationIndex>
523template <
typename Rhs>
525 const Index cols = this->cols();
526 const Index nrhs = rhs.cols();
531 rhs.row(k).swap(rhs.row(
rank - 1));
533 rhs.middleRows(
rank - 1, cols -
rank + 1)
536 rhs.row(k).swap(rhs.row(
rank - 1));
541#ifndef EIGEN_PARSED_BY_DOXYGEN
542template <
typename MatrixType_,
typename PermutationIndex_>
543template <
typename RhsType,
typename DstType>
545 DstType& dst)
const {
546 const Index rank = this->rank();
553 typename RhsType::PlainObject c(rhs);
554 c.applyOnTheLeft(matrixQ().setLength(rank).adjoint());
557 dst.topRows(rank) = matrixT().topLeftCorner(rank, rank).template triangularView<Upper>().solve(c.topRows(rank));
559 const Index cols = this->cols();
563 dst.bottomRows(cols - rank).setZero();
564 applyZAdjointOnTheLeftInPlace(dst);
568 dst = colsPermutation() * dst;
571template <
typename MatrixType_,
typename PermutationIndex_>
572template <
bool Conjugate,
typename RhsType,
typename DstType>
574 DstType& dst)
const {
575 const Index rank = this->rank();
582 typename RhsType::PlainObject c(colsPermutation().transpose() * rhs);
585 applyZOnTheLeftInPlace<!Conjugate>(c);
589 .topLeftCorner(rank, rank)
590 .template triangularView<Upper>()
592 .template conjugateIf<Conjugate>()
593 .solveInPlace(c.topRows(rank));
595 dst.topRows(rank) = c.topRows(rank);
596 dst.bottomRows(rows() - rank).setZero();
598 dst.applyOnTheLeft(householderQ().setLength(rank).
template conjugateIf<!Conjugate>());
604template <
typename MatrixType,
typename PermutationIndex>
605struct traits<Inverse<CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>>>
606 : traits<typename Transpose<typename MatrixType::PlainObject>::PlainObject> {
610template <
typename DstXprType,
typename MatrixType,
typename PermutationIndex>
611struct Assignment<DstXprType, Inverse<CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>>,
612 internal::assign_op<typename DstXprType::Scalar,
613 typename CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>::Scalar>,
615 typedef CompleteOrthogonalDecomposition<MatrixType, PermutationIndex> CodType;
616 typedef Inverse<CodType> SrcXprType;
617 static void run(DstXprType& dst,
const SrcXprType& src,
618 const internal::assign_op<typename DstXprType::Scalar, typename CodType::Scalar>&) {
619 typedef Matrix<
typename CodType::Scalar, CodType::RowsAtCompileTime, CodType::RowsAtCompileTime, 0,
620 CodType::MaxRowsAtCompileTime, CodType::MaxRowsAtCompileTime>
622 dst = src.nestedExpression().solve(IdentityMatrixType::Identity(src.cols(), src.cols()));
629template <
typename MatrixType,
typename PermutationIndex>
630typename CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>::HouseholderSequenceType
632 return m_cpqr.householderQ();
639template <
typename Derived>
640template <
typename PermutationIndex>
642MatrixBase<Derived>::completeOrthogonalDecomposition()
const {
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition ColPivHouseholderQR.h:54
MatrixType::Scalar determinant() const
Definition ColPivHouseholderQR.h:440
Complete orthogonal decomposition (COD) of a matrix.
Definition CompleteOrthogonalDecomposition.h:54
const HCoeffsType & zCoeffs() const
Definition CompleteOrthogonalDecomposition.h:316
MatrixType matrixZ() const
Definition CompleteOrthogonalDecomposition.h:151
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
Definition CompleteOrthogonalDecomposition.h:524
CompleteOrthogonalDecomposition()
Default Constructor.
Definition CompleteOrthogonalDecomposition.h:84
bool isSurjective() const
Definition CompleteOrthogonalDecomposition.h:280
MatrixType::Scalar signDeterminant() const
Definition CompleteOrthogonalDecomposition.h:443
Index dimensionOfKernel() const
Definition CompleteOrthogonalDecomposition.h:262
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition CompleteOrthogonalDecomposition.h:112
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition CompleteOrthogonalDecomposition.h:92
const Solve< CompleteOrthogonalDecomposition, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixType::RealScalar absDeterminant() const
Definition CompleteOrthogonalDecomposition.h:432
HouseholderSequenceType householderQ(void) const
Definition CompleteOrthogonalDecomposition.h:631
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const
Definition CompleteOrthogonalDecomposition.h:296
MatrixType::RealScalar logAbsDeterminant() const
Definition CompleteOrthogonalDecomposition.h:437
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
Definition CompleteOrthogonalDecomposition.h:337
RealScalar threshold() const
Definition CompleteOrthogonalDecomposition.h:359
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was successful.
Definition CompleteOrthogonalDecomposition.h:383
bool isInvertible() const
Definition CompleteOrthogonalDecomposition.h:289
const HCoeffsType & hCoeffs() const
Definition CompleteOrthogonalDecomposition.h:309
const PermutationType & colsPermutation() const
Definition CompleteOrthogonalDecomposition.h:184
const MatrixType & matrixQTZ() const
Definition CompleteOrthogonalDecomposition.h:160
bool isInjective() const
Definition CompleteOrthogonalDecomposition.h:271
void applyZOnTheLeftInPlace(Rhs &rhs) const
Definition CompleteOrthogonalDecomposition.h:504
void computeInPlace()
Definition CompleteOrthogonalDecomposition.h:455
Index nonzeroPivots() const
Definition CompleteOrthogonalDecomposition.h:368
CompleteOrthogonalDecomposition & setThreshold(Default_t)
Definition CompleteOrthogonalDecomposition.h:350
MatrixType::Scalar determinant() const
Definition CompleteOrthogonalDecomposition.h:427
Index rank() const
Definition CompleteOrthogonalDecomposition.h:253
const MatrixType & matrixT() const
Definition CompleteOrthogonalDecomposition.h:173
RealScalar maxPivot() const
Definition CompleteOrthogonalDecomposition.h:373
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition CompleteOrthogonalDecomposition.h:127
EvalReturnType eval() const
Definition DenseBase.h:386
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition HouseholderSequence.h:117
Expression of the inverse of another expression.
Definition Inverse.h:43
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Permutation matrix.
Definition PermutationMatrix.h:283
Pseudo expression representing a solving operation.
Definition Solve.h:62
constexpr Derived & derived()
Definition EigenBase.h:49
SolverBase()
Definition SolverBase.h:105
const ConstTransposeReturnType transpose() const
Definition SolverBase.h:128
const AdjointReturnType adjoint() const
Definition SolverBase.h:144
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Definition EigenBase.h:33
constexpr Derived & derived()
Definition EigenBase.h:49
Eigen::Index Index
The interface type of indices.
Definition EigenBase.h:43