11 #ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H 12 #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H 15 #include "./InternalHeaderCheck.h" 21 template <
typename MatrixType_,
typename PermutationIndex_>
22 struct traits<FullPivHouseholderQR<MatrixType_, PermutationIndex_> > : traits<MatrixType_> {
23 typedef MatrixXpr XprKind;
24 typedef SolverStorage StorageKind;
25 typedef PermutationIndex_ PermutationIndex;
29 template <
typename MatrixType,
typename PermutationIndex>
32 template <
typename MatrixType,
typename PermutationIndex>
34 typedef typename MatrixType::PlainObject ReturnType;
62 template <
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)>
79 IntDiagSizeVectorType;
80 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex> PermutationType;
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;
92 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
104 m_rows_transpositions(),
105 m_cols_transpositions(),
106 m_cols_permutation(),
108 m_isInitialized(false),
109 m_usePrescribedThreshold(false) {}
119 m_hCoeffs((
std::min)(rows, cols)),
120 m_rows_transpositions((
std::min)(rows, cols)),
121 m_cols_transpositions((
std::min)(rows, cols)),
122 m_cols_permutation(cols),
124 m_isInitialized(false),
125 m_usePrescribedThreshold(false) {}
139 template <
typename InputType>
141 : m_qr(matrix.rows(), matrix.cols()),
142 m_hCoeffs((
std::min)(matrix.rows(), matrix.cols())),
143 m_rows_transpositions((
std::min)(matrix.rows(), matrix.cols())),
144 m_cols_transpositions((
std::min)(matrix.rows(), matrix.cols())),
145 m_cols_permutation(matrix.cols()),
146 m_temp(matrix.cols()),
147 m_isInitialized(false),
148 m_usePrescribedThreshold(false) {
159 template <
typename InputType>
162 m_hCoeffs((
std::min)(matrix.rows(), matrix.cols())),
163 m_rows_transpositions((
std::min)(matrix.rows(), matrix.cols())),
164 m_cols_transpositions((
std::min)(matrix.rows(), matrix.cols())),
165 m_cols_permutation(matrix.cols()),
166 m_temp(matrix.cols()),
167 m_isInitialized(false),
168 m_usePrescribedThreshold(false) {
172 #ifdef EIGEN_PARSED_BY_DOXYGEN 188 template <
typename Rhs>
194 MatrixQReturnType
matrixQ(
void)
const;
199 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
203 template <
typename InputType>
208 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
209 return m_cols_permutation;
214 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
215 return m_rows_transpositions;
284 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
285 RealScalar premultiplied_threshold =
abs(m_maxpivot) *
threshold();
287 for (
Index i = 0; i < m_nonzero_pivots; ++i) result += (
abs(m_qr.coeff(i, i)) > premultiplied_threshold);
298 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
299 return cols() -
rank();
310 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
311 return rank() == cols();
322 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
323 return rank() == rows();
333 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
343 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
347 inline Index rows()
const {
return m_qr.rows(); }
348 inline Index cols()
const {
return m_qr.cols(); }
354 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
374 m_usePrescribedThreshold =
true;
388 m_usePrescribedThreshold =
false;
397 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
398 return m_usePrescribedThreshold ? m_prescribedThreshold
412 eigen_assert(m_isInitialized &&
"LU is not initialized.");
413 return m_nonzero_pivots;
421 #ifndef EIGEN_PARSED_BY_DOXYGEN 422 template <
typename RhsType,
typename DstType>
423 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
425 template <
bool Conjugate,
typename RhsType,
typename DstType>
426 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
430 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
432 void computeInPlace();
435 HCoeffsType m_hCoeffs;
436 IntDiagSizeVectorType m_rows_transpositions;
437 IntDiagSizeVectorType m_cols_transpositions;
438 PermutationType m_cols_permutation;
439 RowVectorType m_temp;
440 bool m_isInitialized, m_usePrescribedThreshold;
441 RealScalar m_prescribedThreshold, m_maxpivot;
442 Index m_nonzero_pivots;
443 RealScalar m_precision;
447 template <
typename MatrixType,
typename PermutationIndex>
449 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
450 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
452 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
453 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().prod() : Scalar(0);
456 template <
typename MatrixType,
typename PermutationIndex>
459 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
460 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
461 return isInjective() ?
abs(m_qr.diagonal().prod()) : RealScalar(0);
464 template <
typename MatrixType,
typename PermutationIndex>
466 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
467 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
471 template <
typename MatrixType,
typename PermutationIndex>
473 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
474 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
476 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
477 return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().array().sign().prod() : Scalar(0);
486 template <
typename MatrixType,
typename PermutationIndex>
487 template <
typename InputType>
495 template <
typename MatrixType,
typename PermutationIndex>
499 Index rows = m_qr.rows();
500 Index cols = m_qr.cols();
501 Index size = (std::min)(rows, cols);
503 m_hCoeffs.resize(size);
509 m_rows_transpositions.resize(size);
510 m_cols_transpositions.resize(size);
511 Index number_of_transpositions = 0;
513 RealScalar biggest(0);
515 m_nonzero_pivots = size;
516 m_maxpivot = RealScalar(0);
518 for (
Index k = 0; k < size; ++k) {
519 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
520 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
521 typedef typename Scoring::result_type Score;
523 Score score = m_qr.bottomRightCorner(rows - k, cols - k)
524 .unaryExpr(Scoring())
525 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
526 row_of_biggest_in_corner += k;
527 col_of_biggest_in_corner += k;
528 RealScalar biggest_in_corner =
529 internal::abs_knowing_score<Scalar>()(m_qr(row_of_biggest_in_corner, col_of_biggest_in_corner), score);
530 if (k == 0) biggest = biggest_in_corner;
533 if (internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision)) {
534 m_nonzero_pivots = k;
535 for (
Index i = k; i < size; i++) {
536 m_rows_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
537 m_cols_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
538 m_hCoeffs.coeffRef(i) = Scalar(0);
543 m_rows_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(row_of_biggest_in_corner);
544 m_cols_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(col_of_biggest_in_corner);
545 if (k != row_of_biggest_in_corner) {
546 m_qr.row(k).tail(cols - k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols - k));
547 ++number_of_transpositions;
549 if (k != col_of_biggest_in_corner) {
550 m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
551 ++number_of_transpositions;
555 m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
556 m_qr.coeffRef(k, k) = beta;
559 if (
abs(beta) > m_maxpivot) m_maxpivot =
abs(beta);
561 m_qr.bottomRightCorner(rows - k, cols - k - 1)
562 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows - k - 1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k + 1));
565 m_cols_permutation.setIdentity(cols);
566 for (
Index k = 0; k < size; ++k) m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
568 m_det_p = (number_of_transpositions % 2) ? -1 : 1;
569 m_isInitialized =
true;
572 #ifndef EIGEN_PARSED_BY_DOXYGEN 573 template <
typename MatrixType_,
typename PermutationIndex_>
574 template <
typename RhsType,
typename DstType>
575 void FullPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl(
const RhsType& rhs, DstType& dst)
const {
576 const Index l_rank = rank();
585 typename RhsType::PlainObject c(rhs);
587 Matrix<typename RhsType::Scalar, 1, RhsType::ColsAtCompileTime> temp(rhs.cols());
588 for (
Index k = 0; k < l_rank; ++k) {
589 Index remainingSize = rows() - k;
590 c.row(k).swap(c.row(m_rows_transpositions.coeff(k)));
591 c.bottomRightCorner(remainingSize, rhs.cols())
592 .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
595 m_qr.topLeftCorner(l_rank, l_rank).template triangularView<Upper>().solveInPlace(c.topRows(l_rank));
597 for (
Index i = 0; i < l_rank; ++i) dst.row(m_cols_permutation.indices().coeff(i)) = c.row(i);
598 for (
Index i = l_rank; i < cols(); ++i) dst.row(m_cols_permutation.indices().coeff(i)).setZero();
601 template <
typename MatrixType_,
typename PermutationIndex_>
602 template <
bool Conjugate,
typename RhsType,
typename DstType>
603 void FullPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl_transposed(
const RhsType& rhs,
604 DstType& dst)
const {
605 const Index l_rank = rank();
612 typename RhsType::PlainObject c(m_cols_permutation.transpose() * rhs);
614 m_qr.topLeftCorner(l_rank, l_rank)
615 .template triangularView<Upper>()
617 .template conjugateIf<Conjugate>()
618 .solveInPlace(c.topRows(l_rank));
620 dst.topRows(l_rank) = c.topRows(l_rank);
621 dst.bottomRows(rows() - l_rank).setZero();
623 Matrix<Scalar, 1, DstType::ColsAtCompileTime> temp(dst.cols());
624 const Index size = (std::min)(rows(), cols());
625 for (
Index k = size - 1; k >= 0; --k) {
626 Index remainingSize = rows() - k;
628 dst.bottomRightCorner(remainingSize, dst.cols())
629 .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1).template conjugateIf<!Conjugate>(),
630 m_hCoeffs.template conjugateIf<Conjugate>().coeff(k), &temp.coeffRef(0));
632 dst.row(k).swap(dst.row(m_rows_transpositions.coeff(k)));
639 template <
typename DstXprType,
typename MatrixType,
typename PermutationIndex>
640 struct Assignment<DstXprType, Inverse<FullPivHouseholderQR<MatrixType, PermutationIndex> >,
641 internal::assign_op<typename DstXprType::Scalar,
642 typename FullPivHouseholderQR<MatrixType, PermutationIndex>::Scalar>,
644 typedef FullPivHouseholderQR<MatrixType, PermutationIndex> QrType;
645 typedef Inverse<QrType> SrcXprType;
646 static void run(DstXprType& dst,
const SrcXprType& src,
647 const internal::assign_op<typename DstXprType::Scalar, typename QrType::Scalar>&) {
648 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
658 template <
typename MatrixType,
typename PermutationIndex>
659 struct FullPivHouseholderQRMatrixQReturnType
660 :
public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> > {
662 typedef typename FullPivHouseholderQR<MatrixType, PermutationIndex>::IntDiagSizeVectorType IntDiagSizeVectorType;
663 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
664 typedef Matrix<
typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime,
RowMajor, 1,
665 MatrixType::MaxRowsAtCompileTime>
668 FullPivHouseholderQRMatrixQReturnType(
const MatrixType& qr,
const HCoeffsType& hCoeffs,
669 const IntDiagSizeVectorType& rowsTranspositions)
670 : m_qr(qr), m_hCoeffs(hCoeffs), m_rowsTranspositions(rowsTranspositions) {}
672 template <
typename ResultType>
673 void evalTo(ResultType& result)
const {
674 const Index rows = m_qr.rows();
675 WorkVectorType workspace(rows);
676 evalTo(result, workspace);
679 template <
typename ResultType>
680 void evalTo(ResultType& result, WorkVectorType& workspace)
const {
685 const Index rows = m_qr.rows();
686 const Index cols = m_qr.cols();
687 const Index size = (std::min)(rows, cols);
688 workspace.resize(rows);
689 result.setIdentity(rows, rows);
690 for (
Index k = size - 1; k >= 0; k--) {
691 result.block(k, k, rows - k, rows - k)
692 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows - k - 1),
conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k));
693 result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
697 Index rows()
const {
return m_qr.rows(); }
698 Index cols()
const {
return m_qr.rows(); }
701 typename MatrixType::Nested m_qr;
702 typename HCoeffsType::Nested m_hCoeffs;
703 typename IntDiagSizeVectorType::Nested m_rowsTranspositions;
713 template <
typename MatrixType,
typename PermutationIndex>
714 inline typename FullPivHouseholderQR<MatrixType, PermutationIndex>::MatrixQReturnType
716 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
724 template <
typename Derived>
725 template <
typename PermutationIndex>
733 #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H Index dimensionOfKernel() const
Definition: FullPivHouseholderQR.h:297
RealScalar maxPivot() const
Definition: FullPivHouseholderQR.h:419
bool isSurjective() const
Definition: FullPivHouseholderQR.h:321
Householder rank-revealing QR decomposition of a matrix with full pivoting.
Definition: ForwardDeclarations.h:431
RealScalar threshold() const
Definition: FullPivHouseholderQR.h:396
const MatrixType & matrixQR() const
Definition: FullPivHouseholderQR.h:198
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
MatrixType::Scalar signDeterminant() const
Definition: FullPivHouseholderQR.h:472
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Definition: BFloat16.h:231
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
FullPivHouseholderQR()
Default Constructor.
Definition: FullPivHouseholderQR.h:101
const HCoeffsType & hCoeffs() const
Definition: FullPivHouseholderQR.h:354
Definition: EigenBase.h:33
FullPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: FullPivHouseholderQR.h:373
MatrixQReturnType matrixQ(void) const
Definition: FullPivHouseholderQR.h:715
Expression of the inverse of another expression.
Definition: Inverse.h:43
const IntDiagSizeVectorType & rowsTranspositions() const
Definition: FullPivHouseholderQR.h:213
Index nonzeroPivots() const
Definition: FullPivHouseholderQR.h:411
ComputationInfo info() const
Reports whether the QR factorization was successful.
Definition: FullPivHouseholderQR.h:91
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Expression type for return value of FullPivHouseholderQR::matrixQ()
Definition: FullPivHouseholderQR.h:30
MatrixType::RealScalar logAbsDeterminant() const
Definition: FullPivHouseholderQR.h:465
MatrixType::RealScalar absDeterminant() const
Definition: FullPivHouseholderQR.h:457
const PermutationType & colsPermutation() const
Definition: FullPivHouseholderQR.h:207
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
const Solve< FullPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
FullPivHouseholderQR & setThreshold(Default_t)
Definition: FullPivHouseholderQR.h:387
FullPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: FullPivHouseholderQR.h:117
FullPivHouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: FullPivHouseholderQR.h:140
const Inverse< FullPivHouseholderQR > inverse() const
Definition: FullPivHouseholderQR.h:342
Definition: Constants.h:320
bool isInjective() const
Definition: FullPivHouseholderQR.h:309
Index rank() const
Definition: FullPivHouseholderQR.h:282
Pseudo expression representing a solving operation.
Definition: Solve.h:62
ComputationInfo
Definition: Constants.h:438
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:72
MatrixType::Scalar determinant() const
Definition: FullPivHouseholderQR.h:448
FullPivHouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: FullPivHouseholderQR.h:160
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
bool isInvertible() const
Definition: FullPivHouseholderQR.h:332