10 #ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H 11 #define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H 16 template <
typename _MatrixType>
51 :
public SolverBase<CompleteOrthogonalDecomposition<_MatrixType> >
57 template<
typename Derived>
62 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
63 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
99 : m_cpqr(rows, cols), m_zCoeffs((
std::
min)(rows, cols)), m_temp(cols) {}
117 template <
typename InputType>
119 : m_cpqr(matrix.
rows(), matrix.
cols()),
121 m_temp(matrix.
cols())
132 template<
typename InputType>
134 : m_cpqr(matrix.derived()),
136 m_temp(matrix.
cols())
141 #ifdef EIGEN_PARSED_BY_DOXYGEN 151 template <
typename Rhs>
162 MatrixType
Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
163 applyZOnTheLeftInPlace<false>(
Z);
170 const MatrixType&
matrixQTZ()
const {
return m_cpqr.matrixQR(); }
183 const MatrixType&
matrixT()
const {
return m_cpqr.matrixQR(); }
185 template <
typename InputType>
195 return m_cpqr.colsPermutation();
280 eigen_assert(m_cpqr.m_isInitialized &&
"CompleteOrthogonalDecomposition is not initialized.");
292 inline const HCoeffsType&
hCoeffs()
const {
return m_cpqr.hCoeffs(); }
299 const HCoeffsType&
zCoeffs()
const {
return m_zCoeffs; }
367 eigen_assert(m_cpqr.m_isInitialized &&
"Decomposition is not initialized.");
371 #ifndef EIGEN_PARSED_BY_DOXYGEN 372 template <
typename RhsType,
typename DstType>
373 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
375 template<
bool Conjugate,
typename RhsType,
typename DstType>
376 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
384 template<
bool Transpose_,
typename Rhs>
387 eigen_assert(m_cpqr.m_isInitialized &&
"CompleteOrthogonalDecomposition is not initialized.");
388 eigen_assert((Transpose_?derived().
cols():derived().
rows())==b.rows() &&
"CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
391 void computeInPlace();
397 template <
bool Conjugate,
typename Rhs>
398 void applyZOnTheLeftInPlace(
Rhs& rhs)
const;
402 template <
typename Rhs>
403 void applyZAdjointOnTheLeftInPlace(
Rhs& rhs)
const;
410 template <
typename MatrixType>
413 return m_cpqr.absDeterminant();
416 template <
typename MatrixType>
419 return m_cpqr.logAbsDeterminant();
429 template <
typename MatrixType>
432 check_template_parameters();
437 const Index rank = m_cpqr.rank();
440 m_zCoeffs.resize((
std::min)(rows, cols));
455 for (
Index k = rank - 1; k >= 0; --k) {
460 m_cpqr.m_qr.col(k).head(k + 1).swap(
461 m_cpqr.m_qr.col(rank - 1).head(k + 1));
468 .tail(cols - rank + 1)
469 .makeHouseholderInPlace(m_zCoeffs(k), beta);
470 m_cpqr.m_qr(k, rank - 1) = beta;
473 m_cpqr.m_qr.topRightCorner(k, cols - rank + 1)
474 .applyHouseholderOnTheRight(
475 m_cpqr.m_qr.row(k).tail(cols - rank).adjoint(), m_zCoeffs(k),
480 m_cpqr.m_qr.col(k).head(k + 1).swap(
481 m_cpqr.m_qr.col(rank - 1).head(k + 1));
487 template <
typename MatrixType>
488 template <
bool Conjugate,
typename Rhs>
492 const Index nrhs = rhs.cols();
493 const Index rank = this->rank();
495 for (
Index k = rank-1; k >= 0; --k) {
497 rhs.row(k).swap(rhs.row(rank - 1));
499 rhs.middleRows(rank - 1, cols - rank + 1)
500 .applyHouseholderOnTheLeft(
501 matrixQTZ().
row(k).
tail(cols - rank).transpose().
template conjugateIf<!Conjugate>(), zCoeffs().
template conjugateIf<Conjugate>()(k),
504 rhs.row(k).swap(rhs.row(rank - 1));
509 template <
typename MatrixType>
510 template <
typename Rhs>
514 const Index nrhs = rhs.cols();
515 const Index rank = this->rank();
517 for (
Index k = 0; k < rank; ++k) {
519 rhs.row(k).swap(rhs.row(rank - 1));
521 rhs.middleRows(rank - 1, cols - rank + 1)
522 .applyHouseholderOnTheLeft(
526 rhs.row(k).swap(rhs.row(rank - 1));
531 #ifndef EIGEN_PARSED_BY_DOXYGEN 532 template <
typename _MatrixType>
533 template <
typename RhsType,
typename DstType>
535 const RhsType& rhs, DstType& dst)
const {
536 const Index rank = this->rank();
543 typename RhsType::PlainObject
c(rhs);
544 c.applyOnTheLeft(matrixQ().setLength(rank).
adjoint());
547 dst.topRows(rank) = matrixT()
548 .topLeftCorner(rank, rank)
549 .template triangularView<Upper>()
550 .solve(c.topRows(rank));
556 dst.bottomRows(cols - rank).setZero();
557 applyZAdjointOnTheLeftInPlace(dst);
561 dst = colsPermutation() * dst;
564 template<
typename _MatrixType>
565 template<
bool Conjugate,
typename RhsType,
typename DstType>
568 const Index rank = this->rank();
575 typename RhsType::PlainObject
c(colsPermutation().transpose()*rhs);
578 applyZOnTheLeftInPlace<!Conjugate>(
c);
581 matrixT().topLeftCorner(rank, rank)
582 .template triangularView<Upper>()
583 .transpose().template conjugateIf<Conjugate>()
584 .solveInPlace(c.topRows(rank));
586 dst.topRows(rank) = c.topRows(rank);
587 dst.bottomRows(
rows()-rank).setZero();
589 dst.applyOnTheLeft(householderQ().setLength(rank).
template conjugateIf<!Conjugate>() );
595 template<
typename MatrixType>
597 :
traits<typename Transpose<typename MatrixType::PlainObject>::PlainObject>
602 template<
typename DstXprType,
typename MatrixType>
617 template <
typename MatrixType>
620 return m_cpqr.householderQ();
627 template <
typename Derived>
635 #endif // EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H #define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
const MatrixType & matrixT() const
HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
Inverse< CodType > SrcXprType
CompleteOrthogonalDecomposition()
Default Constructor.
HouseholderSequenceType matrixQ(void) const
void adjoint(const MatrixType &m)
const CompleteOrthogonalDecomposition< PlainObject > completeOrthogonalDecomposition() const
internal::traits< CompleteOrthogonalDecomposition< _MatrixType > >::Scalar Scalar
MatrixType matrixZ() const
MatrixType::PlainObject PlainObject
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Namespace containing all symbols from the Eigen library.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE FixedSegmentReturnType< internal::get_fixed_value< NType >::value >::Type tail(NType n)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
MatrixType::RealScalar logAbsDeterminant() const
bool isSurjective() const
bool isInvertible() const
Complete orthogonal decomposition (COD) of a matrix.
Eigen::Index Index
The interface type of indices.
internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
SolverStorage StorageKind
internal::plain_diag_type< MatrixType >::type HCoeffsType
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
void _solve_impl(const RhsType &rhs, DstType &dst) const
Sequence of Householder reflections acting on subspaces with decreasing size.
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
const PermutationType & colsPermutation() const
Expression of the inverse of another expression.
RealScalar threshold() const
const MatrixType & matrixQTZ() const
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
void applyZOnTheLeftInPlace(Rhs &rhs) const
Index dimensionOfKernel() const
const HCoeffsType & hCoeffs() const
ColPivHouseholderQR< MatrixType > m_cpqr
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
HouseholderSequenceType householderQ(void) const
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationType
NumTraits< Scalar >::Real RealScalar
CompleteOrthogonalDecomposition & compute(const EigenBase< InputType > &matrix)
Index nonzeroPivots() const
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
MatrixType::RealScalar absDeterminant() const
PermutationType::Index PermIndexType
const HCoeffsType & zCoeffs() const
internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
SolverBase< CompleteOrthogonalDecomposition > Base
RealScalar maxPivot() const
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
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)
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const
Generic expression where a coefficient-wise unary operator is applied to an expression.
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.
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was successful.
A base class for matrix decomposition and solvers.
void _check_solve_assertion(const Rhs &b) const
EIGEN_DEVICE_FUNC Derived & derived()
Base class for all dense matrices, vectors, and expressions.
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
EIGEN_DEVICE_FUNC const XprTypeNestedCleaned & nestedExpression() const
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< typename DstXprType::Scalar, typename CodType::Scalar > &)
CompleteOrthogonalDecomposition< MatrixType > CodType
#define EIGEN_ONLY_USED_FOR_DEBUG(x)