11#ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
12#define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
15#include "./InternalHeaderCheck.h"
21template <
typename MatrixType_,
typename PermutationIndex_>
22struct traits<FullPivHouseholderQR<MatrixType_, PermutationIndex_> > : traits<MatrixType_> {
23 typedef MatrixXpr XprKind;
24 typedef SolverStorage StorageKind;
25 typedef PermutationIndex_ PermutationIndex;
29template <
typename MatrixType,
typename PermutationIndex>
30struct FullPivHouseholderQRMatrixQReturnType;
32template <
typename MatrixType,
typename PermutationIndex>
33struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> > {
34 typedef typename MatrixType::PlainObject ReturnType;
62template <
typename MatrixType_,
typename PermutationIndex_>
65 typedef MatrixType_ MatrixType;
68 typedef PermutationIndex_ PermutationIndex;
72 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
75 typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> MatrixQReturnType;
76 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
77 typedef Matrix<PermutationIndex, 1, internal::min_size_prefer_dynamic(ColsAtCompileTime, RowsAtCompileTime),
RowMajor,
78 1, internal::min_size_prefer_fixed(MaxColsAtCompileTime, MaxRowsAtCompileTime)>
81 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
82 typedef typename internal::plain_col_type<MatrixType>::type ColVectorType;
83 typedef typename MatrixType::PlainObject PlainObject;
93 m_rows_transpositions(),
94 m_cols_transpositions(),
97 m_isInitialized(false),
98 m_usePrescribedThreshold(false) {}
108 m_hCoeffs((std::min)(rows, cols)),
109 m_rows_transpositions((std::min)(rows, cols)),
110 m_cols_transpositions((std::min)(rows, cols)),
111 m_cols_permutation(cols),
113 m_isInitialized(false),
114 m_usePrescribedThreshold(false) {}
128 template <
typename InputType>
130 : m_qr(matrix.rows(), matrix.cols()),
131 m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
132 m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
133 m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
134 m_cols_permutation(matrix.cols()),
135 m_temp(matrix.cols()),
136 m_isInitialized(false),
137 m_usePrescribedThreshold(false) {
148 template <
typename InputType>
151 m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
152 m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
153 m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
154 m_cols_permutation(matrix.cols()),
155 m_temp(matrix.cols()),
156 m_isInitialized(false),
157 m_usePrescribedThreshold(false) {
161#ifdef EIGEN_PARSED_BY_DOXYGEN
177 template <
typename Rhs>
183 MatrixQReturnType
matrixQ(
void)
const;
188 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
192 template <
typename InputType>
197 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
198 return m_cols_permutation;
203 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
204 return m_rows_transpositions;
273 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
274 RealScalar premultiplied_threshold =
abs(m_maxpivot) *
threshold();
276 for (
Index i = 0; i < m_nonzero_pivots; ++i) result += (
abs(m_qr.coeff(i, i)) > premultiplied_threshold);
287 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
288 return cols() -
rank();
299 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
300 return rank() == cols();
311 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
312 return rank() == rows();
322 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
332 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
336 inline Index rows()
const {
return m_qr.rows(); }
337 inline Index cols()
const {
return m_qr.cols(); }
343 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
363 m_usePrescribedThreshold =
true;
377 m_usePrescribedThreshold =
false;
386 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
387 return m_usePrescribedThreshold ? m_prescribedThreshold
401 eigen_assert(m_isInitialized &&
"LU is not initialized.");
402 return m_nonzero_pivots;
410#ifndef EIGEN_PARSED_BY_DOXYGEN
411 template <
typename RhsType,
typename DstType>
412 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
414 template <
bool Conjugate,
typename RhsType,
typename DstType>
415 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
419 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
421 void computeInPlace();
424 HCoeffsType m_hCoeffs;
425 IntDiagSizeVectorType m_rows_transpositions;
426 IntDiagSizeVectorType m_cols_transpositions;
427 PermutationType m_cols_permutation;
428 RowVectorType m_temp;
429 bool m_isInitialized, m_usePrescribedThreshold;
430 RealScalar m_prescribedThreshold, m_maxpivot;
431 Index m_nonzero_pivots;
432 RealScalar m_precision;
436template <
typename MatrixType,
typename PermutationIndex>
438 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
439 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
441 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
442 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().prod() : Scalar(0);
445template <
typename MatrixType,
typename PermutationIndex>
448 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
449 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
450 return isInjective() ? abs(m_qr.diagonal().prod()) : RealScalar(0);
453template <
typename MatrixType,
typename PermutationIndex>
455 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
456 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
460template <
typename MatrixType,
typename PermutationIndex>
462 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
463 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
465 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
466 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().array().sign().prod() : Scalar(0);
475template <
typename MatrixType,
typename PermutationIndex>
476template <
typename InputType>
484template <
typename MatrixType,
typename PermutationIndex>
488 Index rows = m_qr.rows();
489 Index cols = m_qr.cols();
490 Index size = (std::min)(rows, cols);
492 m_hCoeffs.resize(size);
498 m_rows_transpositions.resize(size);
499 m_cols_transpositions.resize(size);
500 Index number_of_transpositions = 0;
502 RealScalar biggest(0);
504 m_nonzero_pivots = size;
505 m_maxpivot = RealScalar(0);
507 for (
Index k = 0; k < size; ++k) {
508 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
509 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
510 typedef typename Scoring::result_type Score;
512 Score score = m_qr.bottomRightCorner(rows - k, cols - k)
513 .unaryExpr(Scoring())
514 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
515 row_of_biggest_in_corner += k;
516 col_of_biggest_in_corner += k;
517 RealScalar biggest_in_corner =
518 internal::abs_knowing_score<Scalar>()(m_qr(row_of_biggest_in_corner, col_of_biggest_in_corner), score);
519 if (k == 0) biggest = biggest_in_corner;
522 if (internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision)) {
523 m_nonzero_pivots = k;
524 for (
Index i = k; i < size; i++) {
525 m_rows_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
526 m_cols_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
527 m_hCoeffs.coeffRef(i) = Scalar(0);
532 m_rows_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(row_of_biggest_in_corner);
533 m_cols_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(col_of_biggest_in_corner);
534 if (k != row_of_biggest_in_corner) {
535 m_qr.row(k).tail(cols - k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols - k));
536 ++number_of_transpositions;
538 if (k != col_of_biggest_in_corner) {
539 m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
540 ++number_of_transpositions;
544 m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
545 m_qr.coeffRef(k, k) = beta;
548 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 m_cols_permutation.setIdentity(cols);
555 for (Index k = 0; k < size; ++k) m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
557 m_det_p = (number_of_transpositions % 2) ? -1 : 1;
558 m_isInitialized =
true;
561#ifndef EIGEN_PARSED_BY_DOXYGEN
562template <
typename MatrixType_,
typename PermutationIndex_>
563template <
typename RhsType,
typename DstType>
564void FullPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl(
const RhsType& rhs, DstType& dst)
const {
565 const Index l_rank = rank();
574 typename RhsType::PlainObject c(rhs);
576 Matrix<typename RhsType::Scalar, 1, RhsType::ColsAtCompileTime> temp(rhs.cols());
577 for (Index k = 0; k < l_rank; ++k) {
578 Index remainingSize = rows() - k;
579 c.row(k).swap(c.row(m_rows_transpositions.coeff(k)));
580 c.bottomRightCorner(remainingSize, rhs.cols())
581 .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
584 m_qr.topLeftCorner(l_rank, l_rank).template triangularView<Upper>().solveInPlace(c.topRows(l_rank));
586 for (Index i = 0; i < l_rank; ++i) dst.row(m_cols_permutation.indices().coeff(i)) = c.row(i);
587 for (Index i = l_rank; i < cols(); ++i) dst.row(m_cols_permutation.indices().coeff(i)).setZero();
590template <
typename MatrixType_,
typename PermutationIndex_>
591template <
bool Conjugate,
typename RhsType,
typename DstType>
592void FullPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl_transposed(
const RhsType& rhs,
593 DstType& dst)
const {
594 const Index l_rank = rank();
601 typename RhsType::PlainObject c(m_cols_permutation.transpose() * rhs);
603 m_qr.topLeftCorner(l_rank, l_rank)
604 .template triangularView<Upper>()
606 .template conjugateIf<Conjugate>()
607 .solveInPlace(c.topRows(l_rank));
609 dst.topRows(l_rank) = c.topRows(l_rank);
610 dst.bottomRows(rows() - l_rank).setZero();
612 Matrix<Scalar, 1, DstType::ColsAtCompileTime> temp(dst.cols());
613 const Index size = (std::min)(rows(), cols());
614 for (Index k = size - 1; k >= 0; --k) {
615 Index remainingSize = rows() - k;
617 dst.bottomRightCorner(remainingSize, dst.cols())
618 .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1).template conjugateIf<!Conjugate>(),
619 m_hCoeffs.template conjugateIf<Conjugate>().coeff(k), &temp.coeffRef(0));
621 dst.row(k).swap(dst.row(m_rows_transpositions.coeff(k)));
628template <
typename DstXprType,
typename MatrixType,
typename PermutationIndex>
629struct Assignment<DstXprType, Inverse<FullPivHouseholderQR<MatrixType, PermutationIndex> >,
630 internal::assign_op<typename DstXprType::Scalar,
631 typename FullPivHouseholderQR<MatrixType, PermutationIndex>::Scalar>,
633 typedef FullPivHouseholderQR<MatrixType, PermutationIndex> QrType;
634 typedef Inverse<QrType> SrcXprType;
635 static void run(DstXprType& dst,
const SrcXprType& src,
636 const internal::assign_op<typename DstXprType::Scalar, typename QrType::Scalar>&) {
637 dst = src.nestedExpression().
solve(MatrixType::Identity(src.rows(), src.cols()));
647template <
typename MatrixType,
typename PermutationIndex>
648struct FullPivHouseholderQRMatrixQReturnType
649 :
public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> > {
651 typedef typename FullPivHouseholderQR<MatrixType, PermutationIndex>::IntDiagSizeVectorType IntDiagSizeVectorType;
652 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
653 typedef Matrix<
typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime,
RowMajor, 1,
654 MatrixType::MaxRowsAtCompileTime>
657 FullPivHouseholderQRMatrixQReturnType(
const MatrixType& qr,
const HCoeffsType& hCoeffs,
658 const IntDiagSizeVectorType& rowsTranspositions)
659 : m_qr(qr), m_hCoeffs(hCoeffs), m_rowsTranspositions(rowsTranspositions) {}
661 template <
typename ResultType>
662 void evalTo(ResultType& result)
const {
663 const Index rows = m_qr.rows();
664 WorkVectorType workspace(rows);
665 evalTo(result, workspace);
668 template <
typename ResultType>
669 void evalTo(ResultType& result, WorkVectorType& workspace)
const {
674 const Index rows = m_qr.rows();
675 const Index cols = m_qr.cols();
676 const Index size = (std::min)(rows, cols);
677 workspace.resize(rows);
678 result.setIdentity(rows, rows);
679 for (Index k = size - 1; k >= 0; k--) {
680 result.block(k, k, rows - k, rows - k)
681 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows - k - 1),
conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k));
682 result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
686 Index rows()
const {
return m_qr.rows(); }
687 Index cols()
const {
return m_qr.rows(); }
690 typename MatrixType::Nested m_qr;
691 typename HCoeffsType::Nested m_hCoeffs;
692 typename IntDiagSizeVectorType::Nested m_rowsTranspositions;
702template <
typename MatrixType,
typename PermutationIndex>
703inline typename FullPivHouseholderQR<MatrixType, PermutationIndex>::MatrixQReturnType
705 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
706 return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
713template <
typename Derived>
714template <
typename PermutationIndex>
Householder rank-revealing QR decomposition of a matrix with full pivoting.
Definition FullPivHouseholderQR.h:63
FullPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition FullPivHouseholderQR.h:362
Index dimensionOfKernel() const
Definition FullPivHouseholderQR.h:286
const MatrixType & matrixQR() const
Definition FullPivHouseholderQR.h:187
const Solve< FullPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixQReturnType matrixQ(void) const
Definition FullPivHouseholderQR.h:704
FullPivHouseholderQR & setThreshold(Default_t)
Definition FullPivHouseholderQR.h:376
bool isInjective() const
Definition FullPivHouseholderQR.h:298
bool isInvertible() const
Definition FullPivHouseholderQR.h:321
MatrixType::RealScalar logAbsDeterminant() const
Definition FullPivHouseholderQR.h:454
Index rank() const
Definition FullPivHouseholderQR.h:271
const PermutationType & colsPermutation() const
Definition FullPivHouseholderQR.h:196
bool isSurjective() const
Definition FullPivHouseholderQR.h:310
const IntDiagSizeVectorType & rowsTranspositions() const
Definition FullPivHouseholderQR.h:202
RealScalar maxPivot() const
Definition FullPivHouseholderQR.h:408
MatrixType::Scalar determinant() const
Definition FullPivHouseholderQR.h:437
FullPivHouseholderQR()
Default Constructor.
Definition FullPivHouseholderQR.h:90
const HCoeffsType & hCoeffs() const
Definition FullPivHouseholderQR.h:343
MatrixType::RealScalar absDeterminant() const
Definition FullPivHouseholderQR.h:446
FullPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition FullPivHouseholderQR.h:106
FullPivHouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition FullPivHouseholderQR.h:129
const Inverse< FullPivHouseholderQR > inverse() const
Definition FullPivHouseholderQR.h:331
FullPivHouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition FullPivHouseholderQR.h:149
Index nonzeroPivots() const
Definition FullPivHouseholderQR.h:400
RealScalar threshold() const
Definition FullPivHouseholderQR.h:385
MatrixType::Scalar signDeterminant() const
Definition FullPivHouseholderQR.h:461
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:280
Pseudo expression representing a solving operation.
Definition Solve.h:62
A base class for matrix decomposition and solvers.
Definition SolverBase.h:72
FullPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived()
Definition EigenBase.h:49
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SolverBase.h:106
@ RowMajor
Definition Constants.h:320
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
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:83
Definition EigenBase.h:33
Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523