11#ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
12#define EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
15#include "./InternalHeaderCheck.h"
20template <
typename MatrixType_,
typename PermutationIndex_>
21struct traits<ColPivHouseholderQR<MatrixType_, PermutationIndex_>> : traits<MatrixType_> {
22 typedef MatrixXpr XprKind;
23 typedef SolverStorage StorageKind;
24 typedef PermutationIndex_ PermutationIndex;
53template <
typename MatrixType_,
typename PermutationIndex_>
56 typedef MatrixType_ MatrixType;
59 typedef PermutationIndex_ PermutationIndex;
63 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
66 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
68 typedef typename internal::plain_row_type<MatrixType, PermutationIndex>::type IntRowVectorType;
69 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
70 typedef typename internal::plain_row_type<MatrixType, RealScalar>::type RealRowVectorType;
72 HouseholderSequenceType;
73 typedef typename MatrixType::PlainObject PlainObject;
77 Index diag = numext::mini(rows, cols);
78 m_hCoeffs.resize(diag);
79 m_colsPermutation.resize(cols);
80 m_colsTranspositions.resize(cols);
82 m_colNormsUpdated.resize(cols);
83 m_colNormsDirect.resize(cols);
84 m_isInitialized =
false;
85 m_usePrescribedThreshold =
false;
99 m_colsTranspositions(),
103 m_isInitialized(false),
104 m_usePrescribedThreshold(false) {}
126 template <
typename InputType>
139 template <
typename InputType>
145#ifdef EIGEN_PARSED_BY_DOXYGEN
160 template <
typename Rhs>
165 HouseholderSequenceType matrixQ()
const {
return householderQ(); }
170 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
184 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
188 template <
typename InputType>
193 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
194 return m_colsPermutation;
263 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
264 RealScalar premultiplied_threshold =
abs(m_maxpivot) *
threshold();
266 for (
Index i = 0; i < m_nonzero_pivots; ++i) result += (
abs(m_qr.coeff(i, i)) > premultiplied_threshold);
277 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
278 return cols() -
rank();
289 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
290 return rank() == cols();
301 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
302 return rank() == rows();
312 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
322 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
326 inline Index rows()
const {
return m_qr.rows(); }
327 inline Index cols()
const {
return m_qr.cols(); }
333 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
353 m_usePrescribedThreshold =
true;
367 m_usePrescribedThreshold =
false;
376 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
377 return m_usePrescribedThreshold ? m_prescribedThreshold
380 : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize());
391 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
392 return m_nonzero_pivots;
407 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
411#ifndef EIGEN_PARSED_BY_DOXYGEN
412 template <
typename RhsType,
typename DstType>
413 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
415 template <
bool Conjugate,
typename RhsType,
typename DstType>
416 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
422 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar)
424 void computeInPlace();
427 HCoeffsType m_hCoeffs;
428 PermutationType m_colsPermutation;
429 IntRowVectorType m_colsTranspositions;
430 RowVectorType m_temp;
431 RealRowVectorType m_colNormsUpdated;
432 RealRowVectorType m_colNormsDirect;
433 bool m_isInitialized, m_usePrescribedThreshold;
434 RealScalar m_prescribedThreshold, m_maxpivot;
435 Index m_nonzero_pivots;
439template <
typename MatrixType,
typename PermutationIndex>
441 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
442 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
444 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
445 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().prod() : Scalar(0);
448template <
typename MatrixType,
typename PermutationIndex>
451 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
452 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
453 return isInjective() ?
abs(m_qr.diagonal().prod()) : RealScalar(0);
456template <
typename MatrixType,
typename PermutationIndex>
458 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
459 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
460 return isInjective() ? m_qr.diagonal().cwiseAbs().array().log().sum() : -NumTraits<RealScalar>::infinity();
463template <
typename MatrixType,
typename PermutationIndex>
465 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
466 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
468 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
469 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().array().sign().prod() : Scalar(0);
478template <
typename MatrixType,
typename PermutationIndex>
479template <
typename InputType>
487template <
typename MatrixType,
typename PermutationIndex>
488void ColPivHouseholderQR<MatrixType, PermutationIndex>::computeInPlace() {
489 eigen_assert(m_qr.cols() <= NumTraits<PermutationIndex>::highest());
493 Index rows = m_qr.rows();
494 Index cols = m_qr.cols();
495 Index size = m_qr.diagonalSize();
497 m_hCoeffs.resize(size);
501 m_colsTranspositions.resize(m_qr.cols());
502 Index number_of_transpositions = 0;
504 m_colNormsUpdated.resize(cols);
505 m_colNormsDirect.resize(cols);
506 for (
Index k = 0; k < cols; ++k) {
509 m_colNormsDirect.coeffRef(k) = m_qr.col(k).norm();
510 m_colNormsUpdated.coeffRef(k) = m_colNormsDirect.coeffRef(k);
513 RealScalar threshold_helper =
514 numext::abs2<RealScalar>(m_colNormsUpdated.maxCoeff() * NumTraits<RealScalar>::epsilon()) / RealScalar(rows);
515 RealScalar norm_downdate_threshold = numext::sqrt(NumTraits<RealScalar>::epsilon());
517 m_nonzero_pivots = size;
518 m_maxpivot = RealScalar(0);
520 for (Index k = 0; k < size; ++k) {
522 Index biggest_col_index;
523 RealScalar biggest_col_sq_norm = numext::abs2(m_colNormsUpdated.tail(cols - k).maxCoeff(&biggest_col_index));
524 biggest_col_index += k;
528 if (m_nonzero_pivots == size && biggest_col_sq_norm < threshold_helper * RealScalar(rows - k)) m_nonzero_pivots = k;
531 m_colsTranspositions.coeffRef(k) =
static_cast<PermutationIndex
>(biggest_col_index);
532 if (k != biggest_col_index) {
533 m_qr.col(k).swap(m_qr.col(biggest_col_index));
534 std::swap(m_colNormsUpdated.coeffRef(k), m_colNormsUpdated.coeffRef(biggest_col_index));
535 std::swap(m_colNormsDirect.coeffRef(k), m_colNormsDirect.coeffRef(biggest_col_index));
536 ++number_of_transpositions;
541 m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
544 m_qr.coeffRef(k, k) = beta;
547 if (
abs(beta) > m_maxpivot) m_maxpivot =
abs(beta);
550 m_qr.bottomRightCorner(rows - k, cols - k - 1)
551 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows - k - 1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k + 1));
554 for (
Index j = k + 1; j < cols; ++j) {
559 if (!numext::is_exactly_zero(m_colNormsUpdated.coeffRef(j))) {
560 RealScalar temp =
abs(m_qr.coeffRef(k, j)) / m_colNormsUpdated.coeffRef(j);
561 temp = (RealScalar(1) + temp) * (RealScalar(1) - temp);
562 temp = temp < RealScalar(0) ? RealScalar(0) : temp;
564 temp * numext::abs2<RealScalar>(m_colNormsUpdated.coeffRef(j) / m_colNormsDirect.coeffRef(j));
565 if (temp2 <= norm_downdate_threshold) {
568 m_colNormsDirect.coeffRef(j) = m_qr.col(j).tail(rows - k - 1).norm();
569 m_colNormsUpdated.coeffRef(j) = m_colNormsDirect.coeffRef(j);
571 m_colNormsUpdated.coeffRef(j) *= numext::sqrt(temp);
577 m_colsPermutation.setIdentity(cols);
578 for (
Index k = 0; k < size ; ++k)
579 m_colsPermutation.applyTranspositionOnTheRight(k,
static_cast<Index>(m_colsTranspositions.coeff(k)));
581 m_det_p = (number_of_transpositions % 2) ? -1 : 1;
582 m_isInitialized =
true;
585#ifndef EIGEN_PARSED_BY_DOXYGEN
586template <
typename MatrixType_,
typename PermutationIndex_>
587template <
typename RhsType,
typename DstType>
589 const Index nonzero_pivots = nonzeroPivots();
591 if (nonzero_pivots == 0) {
596 typename RhsType::PlainObject c(rhs);
598 c.applyOnTheLeft(householderQ().setLength(nonzero_pivots).adjoint());
600 m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
601 .template triangularView<Upper>()
602 .solveInPlace(c.topRows(nonzero_pivots));
604 for (
Index i = 0; i < nonzero_pivots; ++i) dst.row(m_colsPermutation.indices().coeff(i)) = c.row(i);
605 for (
Index i = nonzero_pivots; i < cols(); ++i) dst.row(m_colsPermutation.indices().coeff(i)).setZero();
608template <
typename MatrixType_,
typename PermutationIndex_>
609template <
bool Conjugate,
typename RhsType,
typename DstType>
611 DstType& dst)
const {
612 const Index nonzero_pivots = nonzeroPivots();
614 if (nonzero_pivots == 0) {
619 typename RhsType::PlainObject c(m_colsPermutation.transpose() * rhs);
621 m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
622 .template triangularView<Upper>()
624 .template conjugateIf<Conjugate>()
625 .solveInPlace(c.topRows(nonzero_pivots));
627 dst.topRows(nonzero_pivots) = c.topRows(nonzero_pivots);
628 dst.bottomRows(rows() - nonzero_pivots).setZero();
630 dst.applyOnTheLeft(householderQ().setLength(nonzero_pivots).
template conjugateIf<!Conjugate>());
636template <
typename DstXprType,
typename MatrixType,
typename PermutationIndex>
637struct Assignment<DstXprType, Inverse<ColPivHouseholderQR<MatrixType, PermutationIndex>>,
638 internal::assign_op<typename DstXprType::Scalar,
639 typename ColPivHouseholderQR<MatrixType, PermutationIndex>::Scalar>,
641 typedef ColPivHouseholderQR<MatrixType, PermutationIndex> QrType;
642 typedef Inverse<QrType> SrcXprType;
643 static void run(DstXprType& dst,
const SrcXprType& src,
644 const internal::assign_op<typename DstXprType::Scalar, typename QrType::Scalar>&) {
645 dst = src.nestedExpression().
solve(MatrixType::Identity(src.rows(), src.cols()));
654template <
typename MatrixType,
typename PermutationIndex>
655typename ColPivHouseholderQR<MatrixType, PermutationIndex>::HouseholderSequenceType
657 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
658 return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
665template <
typename Derived>
666template <
typename PermutationIndexType>
668MatrixBase<Derived>::colPivHouseholderQr()
const {
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition ColPivHouseholderQR.h:54
RealScalar maxPivot() const
Definition ColPivHouseholderQR.h:398
const Inverse< ColPivHouseholderQR > inverse() const
Definition ColPivHouseholderQR.h:321
const HCoeffsType & hCoeffs() const
Definition ColPivHouseholderQR.h:333
MatrixType::Scalar signDeterminant() const
Definition ColPivHouseholderQR.h:464
Index rank() const
Definition ColPivHouseholderQR.h:261
const Solve< ColPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixType::RealScalar absDeterminant() const
Definition ColPivHouseholderQR.h:449
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition ColPivHouseholderQR.h:352
bool isSurjective() const
Definition ColPivHouseholderQR.h:300
Index nonzeroPivots() const
Definition ColPivHouseholderQR.h:390
bool isInjective() const
Definition ColPivHouseholderQR.h:288
ColPivHouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition ColPivHouseholderQR.h:140
ColPivHouseholderQR & setThreshold(Default_t)
Definition ColPivHouseholderQR.h:366
bool isInvertible() const
Definition ColPivHouseholderQR.h:311
const PermutationType & colsPermutation() const
Definition ColPivHouseholderQR.h:192
const MatrixType & matrixQR() const
Definition ColPivHouseholderQR.h:169
Index dimensionOfKernel() const
Definition ColPivHouseholderQR.h:276
ComputationInfo info() const
Reports whether the QR factorization was successful.
Definition ColPivHouseholderQR.h:406
ColPivHouseholderQR()
Default Constructor.
Definition ColPivHouseholderQR.h:95
HouseholderSequenceType householderQ() const
Definition ColPivHouseholderQR.h:656
MatrixType::Scalar determinant() const
Definition ColPivHouseholderQR.h:440
MatrixType::RealScalar logAbsDeterminant() const
Definition ColPivHouseholderQR.h:457
RealScalar threshold() const
Definition ColPivHouseholderQR.h:375
ColPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition ColPivHouseholderQR.h:112
const MatrixType & matrixR() const
Definition ColPivHouseholderQR.h:183
ColPivHouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition ColPivHouseholderQR.h:127
Complete orthogonal decomposition (COD) of a matrix.
Definition CompleteOrthogonalDecomposition.h:54
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
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 Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SolverBase.h:114
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Definition EigenBase.h:33
constexpr Index cols() const noexcept
Definition EigenBase.h:61
constexpr Derived & derived()
Definition EigenBase.h:49
Eigen::Index Index
The interface type of indices.
Definition EigenBase.h:43
constexpr Index rows() const noexcept
Definition EigenBase.h:59