10#ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
11#define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
16template <
typename _MatrixType>
17struct traits<CompleteOrthogonalDecomposition<_MatrixType> >
18 : traits<_MatrixType> {
47template <
typename _MatrixType>
50 typedef _MatrixType MatrixType;
52 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
53 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
54 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
55 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
57 typedef typename MatrixType::Scalar Scalar;
58 typedef typename MatrixType::RealScalar RealScalar;
59 typedef typename MatrixType::StorageIndex StorageIndex;
60 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
63 typedef typename internal::plain_row_type<MatrixType, Index>::type
65 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
66 typedef typename internal::plain_row_type<MatrixType, RealScalar>::type
69 MatrixType,
typename internal::remove_all<
70 typename HCoeffsType::ConjugateReturnType>::type>
71 HouseholderSequenceType;
72 typedef typename MatrixType::PlainObject PlainObject;
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()),
115 m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
116 m_temp(matrix.cols())
127 template<
typename InputType>
129 : m_cpqr(matrix.derived()),
130 m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
131 m_temp(matrix.cols())
146 template <
typename Rhs>
149 eigen_assert(m_cpqr.m_isInitialized &&
150 "CompleteOrthogonalDecomposition is not initialized.");
155 HouseholderSequenceType matrixQ(
void)
const {
return m_cpqr.householderQ(); }
160 MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
168 const MatrixType&
matrixQTZ()
const {
return m_cpqr.matrixQR(); }
181 const MatrixType&
matrixT()
const {
return m_cpqr.matrixQR(); }
183 template <
typename InputType>
186 m_cpqr.compute(matrix);
193 return m_cpqr.colsPermutation();
281 inline Index rows()
const {
return m_cpqr.rows(); }
282 inline Index cols()
const {
return m_cpqr.cols(); }
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;
374 static void check_template_parameters() {
375 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar);
382 template <
typename Rhs>
386 HCoeffsType m_zCoeffs;
387 RowVectorType m_temp;
390template <
typename MatrixType>
391typename MatrixType::RealScalar
396template <
typename MatrixType>
397typename MatrixType::RealScalar
399 return m_cpqr.logAbsDeterminant();
409template <
typename MatrixType>
412 check_template_parameters();
415 eigen_assert(m_cpqr.cols() <= NumTraits<int>::highest());
418 const Index cols = m_cpqr.cols();
419 const Index rows = m_cpqr.rows();
420 m_zCoeffs.resize((std::min)(rows, cols));
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));
467template <
typename MatrixType>
468template <
typename Rhs>
471 const Index cols = this->cols();
472 const Index nrhs = rhs.cols();
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
490template <
typename _MatrixType>
491template <
typename RhsType,
typename DstType>
493 const RhsType& rhs, DstType& dst)
const {
494 eigen_assert(rhs.rows() == this->rows());
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));
515 const Index cols = this->cols();
519 dst.bottomRows(cols - rank).setZero();
520 applyZAdjointOnTheLeftInPlace(dst);
524 dst = colsPermutation() * dst;
530template<
typename DstXprType,
typename MatrixType>
531struct Assignment<DstXprType, Inverse<CompleteOrthogonalDecomposition<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename CompleteOrthogonalDecomposition<MatrixType>::Scalar>, Dense2Dense>
533 typedef CompleteOrthogonalDecomposition<MatrixType> CodType;
534 typedef Inverse<CodType> SrcXprType;
535 static EIGEN_DEVICE_FUNC
void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename CodType::Scalar> &)
537 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.rows()));
544template <
typename MatrixType>
545typename CompleteOrthogonalDecomposition<MatrixType>::HouseholderSequenceType
547 return m_cpqr.householderQ();
554template <
typename Derived>
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition ColPivHouseholderQR.h:49
MatrixType::RealScalar absDeterminant() const
Definition ColPivHouseholderQR.h:448
Complete orthogonal decomposition (COD) of a matrix.
Definition CompleteOrthogonalDecomposition.h:48
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition CompleteOrthogonalDecomposition.h:128
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
Definition CompleteOrthogonalDecomposition.h:469
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const
Definition CompleteOrthogonalDecomposition.h:276
const HCoeffsType & zCoeffs() const
Definition CompleteOrthogonalDecomposition.h:296
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was succesful.
Definition CompleteOrthogonalDecomposition.h:363
bool isInjective() const
Definition CompleteOrthogonalDecomposition.h:251
RealScalar threshold() const
Definition CompleteOrthogonalDecomposition.h:339
CompleteOrthogonalDecomposition & setThreshold(Default_t)
Definition CompleteOrthogonalDecomposition.h:330
MatrixType matrixZ() const
Definition CompleteOrthogonalDecomposition.h:159
bool isSurjective() const
Definition CompleteOrthogonalDecomposition.h:260
RealScalar maxPivot() const
Definition CompleteOrthogonalDecomposition.h:353
CompleteOrthogonalDecomposition()
Default Constructor.
Definition CompleteOrthogonalDecomposition.h:85
bool isInvertible() const
Definition CompleteOrthogonalDecomposition.h:269
const MatrixType & matrixQTZ() const
Definition CompleteOrthogonalDecomposition.h:168
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition CompleteOrthogonalDecomposition.h:93
const PermutationType & colsPermutation() const
Definition CompleteOrthogonalDecomposition.h:192
const MatrixType & matrixT() const
Definition CompleteOrthogonalDecomposition.h:181
MatrixType::RealScalar absDeterminant() const
Definition CompleteOrthogonalDecomposition.h:392
const HCoeffsType & hCoeffs() const
Definition CompleteOrthogonalDecomposition.h:289
HouseholderSequenceType householderQ(void) const
Definition CompleteOrthogonalDecomposition.h:546
Index dimensionOfKernel() const
Definition CompleteOrthogonalDecomposition.h:242
MatrixType::RealScalar logAbsDeterminant() const
Definition CompleteOrthogonalDecomposition.h:398
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
Definition CompleteOrthogonalDecomposition.h:317
void computeInPlace()
Definition CompleteOrthogonalDecomposition.h:410
Index rank() const
Definition CompleteOrthogonalDecomposition.h:233
Index nonzeroPivots() const
Definition CompleteOrthogonalDecomposition.h:348
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition CompleteOrthogonalDecomposition.h:113
const Solve< CompleteOrthogonalDecomposition, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition CompleteOrthogonalDecomposition.h:147
EvalReturnType eval() const
Definition DenseBase.h:401
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition HouseholderSequence.h:121
Expression of the inverse of another expression.
Definition Inverse.h:44
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
const CompleteOrthogonalDecomposition< PlainObject > completeOrthogonalDecomposition() const
Definition CompleteOrthogonalDecomposition.h:556
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
Permutation matrix.
Definition PermutationMatrix.h:298
Pseudo expression representing a solving operation.
Definition Solve.h:63
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
Definition HouseholderSequence.h:452
ComputationInfo
Definition Constants.h:430
@ Success
Definition Constants.h:432
Namespace containing all symbols from the Eigen library.
Definition A05_PortingFrom2To3.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:65
Definition EigenBase.h:30
Eigen::Index Index
Definition EigenBase.h:38
Derived & derived()
Definition EigenBase.h:45