10 #ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
11 #define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
18 template <
typename MatrixType_,
typename PermutationIndex_>
19 struct traits<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_> >
20 : traits<MatrixType_> {
21 typedef MatrixXpr XprKind;
22 typedef SolverStorage StorageKind;
23 typedef PermutationIndex_ PermutationIndex;
53 :
public SolverBase<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_> >
59 template<
typename Derived>
60 friend struct internal::solve_assertion;
67 typedef typename internal::plain_diag_type<MatrixType>::type
HCoeffsType;
70 typedef typename internal::plain_row_type<MatrixType, Index>::type
72 typedef typename internal::plain_row_type<MatrixType>::type
RowVectorType;
73 typedef typename internal::plain_row_type<MatrixType, RealScalar>::type
77 typename HCoeffsType::ConjugateReturnType>>
116 template <
typename InputType>
131 template<
typename InputType>
140 #ifdef EIGEN_PARSED_BY_DOXYGEN
150 template <
typename Rhs>
162 applyZOnTheLeftInPlace<false>(Z);
184 template <
typename InputType>
385 #ifndef EIGEN_PARSED_BY_DOXYGEN
386 template <
typename RhsType,
typename DstType>
387 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
389 template<
bool Conjugate,
typename RhsType,
typename DstType>
390 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
396 template<
bool Transpose_,
typename Rhs>
400 eigen_assert((Transpose_?
derived().
cols():
derived().
rows())==
b.rows() &&
"CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
409 template <
bool Conjugate,
typename Rhs>
414 template <
typename Rhs>
422 template <
typename MatrixType,
typename PermutationIndex>
425 return m_cpqr.determinant();
428 template <
typename MatrixType,
typename PermutationIndex>
431 return m_cpqr.absDeterminant();
434 template <
typename MatrixType,
typename PermutationIndex>
437 return m_cpqr.logAbsDeterminant();
447 template <
typename MatrixType,
typename PermutationIndex>
452 const Index rank = m_cpqr.rank();
470 for (
Index k = rank - 1; k >= 0; --k) {
475 m_cpqr.m_qr.col(k).head(k + 1).swap(
476 m_cpqr.m_qr.col(rank - 1).head(k + 1));
483 .tail(
cols - rank + 1)
484 .makeHouseholderInPlace(m_zCoeffs(k), beta);
485 m_cpqr.m_qr(k, rank - 1) = beta;
488 m_cpqr.m_qr.topRightCorner(k,
cols - rank + 1)
489 .applyHouseholderOnTheRight(
490 m_cpqr.m_qr.row(k).tail(
cols - rank).adjoint(), m_zCoeffs(k),
495 m_cpqr.m_qr.col(k).head(k + 1).swap(
496 m_cpqr.m_qr.col(rank - 1).head(k + 1));
502 template <
typename MatrixType,
typename PermutationIndex>
503 template <
bool Conjugate,
typename Rhs>
507 const Index nrhs = rhs.cols();
508 const Index rank = this->rank();
510 for (
Index k = rank-1; k >= 0; --k) {
512 rhs.row(k).swap(rhs.row(rank - 1));
514 rhs.middleRows(rank - 1,
cols - rank + 1)
515 .applyHouseholderOnTheLeft(
516 matrixQTZ().
row(k).
tail(
cols - rank).transpose().
template conjugateIf<!Conjugate>(), zCoeffs().
template conjugateIf<Conjugate>()(k),
519 rhs.row(k).swap(rhs.row(rank - 1));
524 template <
typename MatrixType,
typename PermutationIndex>
525 template <
typename Rhs>
529 const Index nrhs = rhs.cols();
530 const Index rank = this->rank();
532 for (
Index k = 0; k < rank; ++k) {
534 rhs.row(k).swap(rhs.row(rank - 1));
536 rhs.middleRows(rank - 1,
cols - rank + 1)
537 .applyHouseholderOnTheLeft(
538 matrixQTZ().
row(k).
tail(
cols - rank).adjoint(), zCoeffs()(k),
541 rhs.row(k).swap(rhs.row(rank - 1));
546 #ifndef EIGEN_PARSED_BY_DOXYGEN
547 template <
typename MatrixType_,
typename PermutationIndex_>
548 template <
typename RhsType,
typename DstType>
550 const RhsType& rhs, DstType& dst)
const {
551 const Index rank = this->rank();
558 typename RhsType::PlainObject
c(rhs);
559 c.applyOnTheLeft(matrixQ().setLength(rank).adjoint());
562 dst.topRows(rank) = matrixT()
563 .topLeftCorner(rank, rank)
564 .template triangularView<Upper>()
565 .solve(
c.topRows(rank));
571 dst.bottomRows(
cols - rank).setZero();
572 applyZAdjointOnTheLeftInPlace(dst);
576 dst = colsPermutation() * dst;
579 template<
typename MatrixType_,
typename PermutationIndex_>
580 template<
bool Conjugate,
typename RhsType,
typename DstType>
581 void CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const
583 const Index rank = this->rank();
590 typename RhsType::PlainObject
c(colsPermutation().transpose()*rhs);
593 applyZOnTheLeftInPlace<!Conjugate>(
c);
596 matrixT().topLeftCorner(rank, rank)
597 .template triangularView<Upper>()
598 .transpose().template conjugateIf<Conjugate>()
599 .solveInPlace(
c.topRows(rank));
601 dst.topRows(rank) =
c.topRows(rank);
604 dst.applyOnTheLeft(householderQ().setLength(rank).
template conjugateIf<!Conjugate>() );
610 template<
typename MatrixType,
typename PermutationIndex>
611 struct traits<Inverse<CompleteOrthogonalDecomposition<
MatrixType, PermutationIndex> > >
612 : traits<typename Transpose<typename MatrixType::PlainObject>::PlainObject>
617 template<
typename DstXprType,
typename MatrixType,
typename PermutationIndex>
618 struct Assignment<DstXprType, Inverse<CompleteOrthogonalDecomposition<
MatrixType, PermutationIndex> >,
internal::assign_op<typename DstXprType::Scalar,typename CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>::Scalar>, Dense2Dense>
620 typedef CompleteOrthogonalDecomposition<MatrixType, PermutationIndex> CodType;
621 typedef Inverse<CodType> SrcXprType;
622 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename CodType::Scalar> &)
624 typedef Matrix<typename CodType::Scalar, CodType::RowsAtCompileTime, CodType::RowsAtCompileTime, 0, CodType::MaxRowsAtCompileTime, CodType::MaxRowsAtCompileTime> IdentityMatrixType;
625 dst = src.nestedExpression().solve(IdentityMatrixType::Identity(src.cols(), src.cols()));
632 template <
typename MatrixType,
typename PermutationIndex>
635 return m_cpqr.householderQ();
642 template <
typename Derived>
643 template <
typename PermutationIndex>
RowXpr row(Index i)
This is the const version of row(). */.
FixedSegmentReturnType<... >::Type tail(NType n)
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Matrix< float, 1, Dynamic > MatrixType
RealScalar maxPivot() const
bool isSurjective() const
Index nonzeroPivots() const
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)
const HCoeffsType & hCoeffs() const
bool isInvertible() const
Index dimensionOfKernel() const
const PermutationType & colsPermutation() const
HouseholderSequenceType householderQ() const
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)
RealScalar threshold() const
const MatrixType & matrixQR() const
Complete orthogonal decomposition (COD) of a matrix.
internal::plain_diag_type< MatrixType >::type HCoeffsType
MatrixType matrixZ() const
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
CompleteOrthogonalDecomposition()
Default Constructor.
bool isSurjective() const
Index dimensionOfKernel() const
HouseholderSequenceType matrixQ(void) const
CompleteOrthogonalDecomposition & compute(const EigenBase< InputType > &matrix)
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex > PermutationType
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
const Solve< CompleteOrthogonalDecomposition, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixType::RealScalar absDeterminant() const
const HCoeffsType & zCoeffs() const
HouseholderSequenceType householderQ(void) const
MatrixType::RealScalar logAbsDeterminant() const
RealScalar threshold() const
ColPivHouseholderQR< MatrixType, PermutationIndex > m_cpqr
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was successful.
const PermutationType & colsPermutation() const
internal::plain_row_type< MatrixType >::type RowVectorType
SolverBase< CompleteOrthogonalDecomposition > Base
internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
bool isInvertible() const
const HCoeffsType & hCoeffs() const
void _check_solve_assertion(const Rhs &b) const
MatrixType::PlainObject PlainObject
internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
CompleteOrthogonalDecomposition & setThreshold(Default_t)
void applyZOnTheLeftInPlace(Rhs &rhs) const
const MatrixType & matrixQTZ() const
Index nonzeroPivots() const
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const
MatrixType::Scalar determinant() const
RealScalar maxPivot() const
const MatrixType & matrixT() const
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
PermutationIndex_ PermutationIndex
NumTraits< Scalar >::Real RealScalar
internal::traits< Derived >::Scalar Scalar
Sequence of Householder reflections acting on subspaces with decreasing size.
Expression of the inverse of another expression.
Base class for all dense matrices, vectors, and expressions.
const CompleteOrthogonalDecomposition< PlainObject, PermutationIndex > completeOrthogonalDecomposition() const
The matrix class, also used for vectors and row-vectors.
Base::PlainObject PlainObject
Derived & setZero(Index size)
Pseudo expression representing a solving operation.
A base class for matrix decomposition and solvers.
internal::traits< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > >::Scalar Scalar
CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > & derived()
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
typename remove_all< T >::type remove_all_t
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Eigen::Index Index
The interface type of indices.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.