$darkmode
Eigen  5.0.1-dev
ColPivHouseholderQR.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
12 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 template <typename MatrixType_, typename PermutationIndex_>
21 struct traits<ColPivHouseholderQR<MatrixType_, PermutationIndex_>> : traits<MatrixType_> {
22  typedef MatrixXpr XprKind;
23  typedef SolverStorage StorageKind;
24  typedef PermutationIndex_ PermutationIndex;
25  enum { Flags = 0 };
26 };
27 
28 } // end namespace internal
29 
53 template <typename MatrixType_, typename PermutationIndex_>
54 class ColPivHouseholderQR : public SolverBase<ColPivHouseholderQR<MatrixType_, PermutationIndex_>> {
55  public:
56  typedef MatrixType_ MatrixType;
57  typedef SolverBase<ColPivHouseholderQR> Base;
58  friend class SolverBase<ColPivHouseholderQR>;
59  typedef PermutationIndex_ PermutationIndex;
60  EIGEN_GENERIC_PUBLIC_INTERFACE(ColPivHouseholderQR)
61 
62  enum {
63  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
65  };
66  typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
67  typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex> PermutationType;
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;
71  typedef HouseholderSequence<MatrixType, internal::remove_all_t<typename HCoeffsType::ConjugateReturnType>>
72  HouseholderSequenceType;
73  typedef typename MatrixType::PlainObject PlainObject;
74 
75  private:
76  void init(Index rows, Index cols) {
77  Index diag = numext::mini(rows, cols);
78  m_hCoeffs.resize(diag);
79  m_colsPermutation.resize(cols);
80  m_colsTranspositions.resize(cols);
81  m_temp.resize(cols);
82  m_colNormsUpdated.resize(cols);
83  m_colNormsDirect.resize(cols);
84  m_isInitialized = false;
85  m_usePrescribedThreshold = false;
86  }
87 
88  public:
96  : m_qr(),
97  m_hCoeffs(),
98  m_colsPermutation(),
99  m_colsTranspositions(),
100  m_temp(),
101  m_colNormsUpdated(),
102  m_colNormsDirect(),
103  m_isInitialized(false),
104  m_usePrescribedThreshold(false) {}
105 
112  ColPivHouseholderQR(Index rows, Index cols) : m_qr(rows, cols) { init(rows, cols); }
113 
126  template <typename InputType>
127  explicit ColPivHouseholderQR(const EigenBase<InputType>& matrix) : m_qr(matrix.rows(), matrix.cols()) {
128  init(matrix.rows(), matrix.cols());
129  compute(matrix.derived());
130  }
131 
139  template <typename InputType>
140  explicit ColPivHouseholderQR(EigenBase<InputType>& matrix) : m_qr(matrix.derived()) {
141  init(matrix.rows(), matrix.cols());
142  computeInPlace();
143  }
144 
145 #ifdef EIGEN_PARSED_BY_DOXYGEN
146 
160  template <typename Rhs>
161  inline const Solve<ColPivHouseholderQR, Rhs> solve(const MatrixBase<Rhs>& b) const;
162 #endif
163 
164  HouseholderSequenceType householderQ() const;
165  HouseholderSequenceType matrixQ() const { return householderQ(); }
166 
169  const MatrixType& matrixQR() const {
170  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
171  return m_qr;
172  }
173 
183  const MatrixType& matrixR() const {
184  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
185  return m_qr;
186  }
187 
188  template <typename InputType>
189  ColPivHouseholderQR& compute(const EigenBase<InputType>& matrix);
190 
193  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
194  return m_colsPermutation;
195  }
196 
210  typename MatrixType::Scalar determinant() const;
211 
225  typename MatrixType::RealScalar absDeterminant() const;
226 
239  typename MatrixType::RealScalar logAbsDeterminant() const;
240 
253  typename MatrixType::Scalar signDeterminant() const;
254 
261  inline Index rank() const {
262  using std::abs;
263  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
264  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
265  Index result = 0;
266  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_qr.coeff(i, i)) > premultiplied_threshold);
267  return result;
268  }
269 
276  inline Index dimensionOfKernel() const {
277  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
278  return cols() - rank();
279  }
280 
288  inline bool isInjective() const {
289  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
290  return rank() == cols();
291  }
292 
300  inline bool isSurjective() const {
301  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
302  return rank() == rows();
303  }
304 
311  inline bool isInvertible() const {
312  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
313  return isInjective() && isSurjective();
314  }
315 
321  inline const Inverse<ColPivHouseholderQR> inverse() const {
322  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
323  return Inverse<ColPivHouseholderQR>(*this);
324  }
325 
326  inline Index rows() const { return m_qr.rows(); }
327  inline Index cols() const { return m_qr.cols(); }
328 
333  const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
334 
353  m_usePrescribedThreshold = true;
354  m_prescribedThreshold = threshold;
355  return *this;
356  }
357 
367  m_usePrescribedThreshold = false;
368  return *this;
369  }
370 
375  RealScalar threshold() const {
376  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
377  return m_usePrescribedThreshold ? m_prescribedThreshold
378  // this formula comes from experimenting (see "LU precision tuning" thread on the
379  // list) and turns out to be identical to Higham's formula used already in LDLt.
380  : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize());
381  }
382 
390  inline Index nonzeroPivots() const {
391  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
392  return m_nonzero_pivots;
393  }
394 
398  RealScalar maxPivot() const { return m_maxpivot; }
399 
407  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
408  return Success;
409  }
410 
411 #ifndef EIGEN_PARSED_BY_DOXYGEN
412  template <typename RhsType, typename DstType>
413  void _solve_impl(const RhsType& rhs, DstType& dst) const;
414 
415  template <bool Conjugate, typename RhsType, typename DstType>
416  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
417 #endif
418 
419  protected:
420  friend class CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>;
421 
422  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
423 
424  void computeInPlace();
425 
426  MatrixType m_qr;
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;
436  Index m_det_p;
437 };
438 
439 template <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!");
443  Scalar detQ;
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);
446 }
447 
448 template <typename MatrixType, typename PermutationIndex>
450  using std::abs;
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);
454 }
455 
456 template <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();
461 }
462 
463 template <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!");
467  Scalar detQ;
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);
470 }
471 
478 template <typename MatrixType, typename PermutationIndex>
479 template <typename InputType>
481  const EigenBase<InputType>& matrix) {
482  m_qr = matrix.derived();
483  computeInPlace();
484  return *this;
485 }
486 
487 template <typename MatrixType, typename PermutationIndex>
489  eigen_assert(m_qr.cols() <= NumTraits<PermutationIndex>::highest());
490 
491  using std::abs;
492 
493  Index rows = m_qr.rows();
494  Index cols = m_qr.cols();
495  Index size = m_qr.diagonalSize();
496 
497  m_hCoeffs.resize(size);
498 
499  m_temp.resize(cols);
500 
501  m_colsTranspositions.resize(m_qr.cols());
502  Index number_of_transpositions = 0;
503 
504  m_colNormsUpdated.resize(cols);
505  m_colNormsDirect.resize(cols);
506  for (Index k = 0; k < cols; ++k) {
507  // colNormsDirect(k) caches the most recent directly computed norm of
508  // column k.
509  m_colNormsDirect.coeffRef(k) = m_qr.col(k).norm();
510  m_colNormsUpdated.coeffRef(k) = m_colNormsDirect.coeffRef(k);
511  }
512 
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());
516 
517  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
518  m_maxpivot = RealScalar(0);
519 
520  for (Index k = 0; k < size; ++k) {
521  // first, we look up in our table m_colNormsUpdated which column has the biggest norm
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;
525 
526  // Track the number of meaningful pivots but do not stop the decomposition to make
527  // sure that the initial matrix is properly reproduced. See bug 941.
528  if (m_nonzero_pivots == size && biggest_col_sq_norm < threshold_helper * RealScalar(rows - k)) m_nonzero_pivots = k;
529 
530  // apply the transposition to the columns
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;
537  }
538 
539  // generate the householder vector, store it below the diagonal
540  RealScalar beta;
541  m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
542 
543  // apply the householder transformation to the diagonal coefficient
544  m_qr.coeffRef(k, k) = beta;
545 
546  // remember the maximum absolute value of diagonal coefficients
547  if (abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
548 
549  // apply the householder transformation
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));
552 
553  // update our table of norms of the columns
554  for (Index j = k + 1; j < cols; ++j) {
555  // The following implements the stable norm downgrade step discussed in
556  // http://www.netlib.org/lapack/lawnspdf/lawn176.pdf
557  // and used in LAPACK routines xGEQPF and xGEQP3.
558  // See lines 278-297 in http://www.netlib.org/lapack/explore-html/dc/df4/sgeqpf_8f_source.html
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;
563  RealScalar temp2 =
564  temp * numext::abs2<RealScalar>(m_colNormsUpdated.coeffRef(j) / m_colNormsDirect.coeffRef(j));
565  if (temp2 <= norm_downdate_threshold) {
566  // The updated norm has become too inaccurate so re-compute the column
567  // norm directly.
568  m_colNormsDirect.coeffRef(j) = m_qr.col(j).tail(rows - k - 1).norm();
569  m_colNormsUpdated.coeffRef(j) = m_colNormsDirect.coeffRef(j);
570  } else {
571  m_colNormsUpdated.coeffRef(j) *= numext::sqrt(temp);
572  }
573  }
574  }
575  }
576 
577  m_colsPermutation.setIdentity(cols);
578  for (Index k = 0; k < size /*m_nonzero_pivots*/; ++k)
579  m_colsPermutation.applyTranspositionOnTheRight(k, static_cast<Index>(m_colsTranspositions.coeff(k)));
580 
581  m_det_p = (number_of_transpositions % 2) ? -1 : 1;
582  m_isInitialized = true;
583 }
584 
585 #ifndef EIGEN_PARSED_BY_DOXYGEN
586 template <typename MatrixType_, typename PermutationIndex_>
587 template <typename RhsType, typename DstType>
588 void ColPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl(const RhsType& rhs, DstType& dst) const {
589  const Index nonzero_pivots = nonzeroPivots();
590 
591  if (nonzero_pivots == 0) {
592  dst.setZero();
593  return;
594  }
595 
596  typename RhsType::PlainObject c(rhs);
597 
598  c.applyOnTheLeft(householderQ().setLength(nonzero_pivots).adjoint());
599 
600  m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
601  .template triangularView<Upper>()
602  .solveInPlace(c.topRows(nonzero_pivots));
603 
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();
606 }
607 
608 template <typename MatrixType_, typename PermutationIndex_>
609 template <bool Conjugate, typename RhsType, typename DstType>
610 void ColPivHouseholderQR<MatrixType_, PermutationIndex_>::_solve_impl_transposed(const RhsType& rhs,
611  DstType& dst) const {
612  const Index nonzero_pivots = nonzeroPivots();
613 
614  if (nonzero_pivots == 0) {
615  dst.setZero();
616  return;
617  }
618 
619  typename RhsType::PlainObject c(m_colsPermutation.transpose() * rhs);
620 
621  m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
622  .template triangularView<Upper>()
623  .transpose()
624  .template conjugateIf<Conjugate>()
625  .solveInPlace(c.topRows(nonzero_pivots));
626 
627  dst.topRows(nonzero_pivots) = c.topRows(nonzero_pivots);
628  dst.bottomRows(rows() - nonzero_pivots).setZero();
629 
630  dst.applyOnTheLeft(householderQ().setLength(nonzero_pivots).template conjugateIf<!Conjugate>());
631 }
632 #endif
633 
634 namespace internal {
635 
636 template <typename DstXprType, typename MatrixType, typename PermutationIndex>
637 struct Assignment<DstXprType, Inverse<ColPivHouseholderQR<MatrixType, PermutationIndex>>,
638  internal::assign_op<typename DstXprType::Scalar,
639  typename ColPivHouseholderQR<MatrixType, PermutationIndex>::Scalar>,
640  Dense2Dense> {
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()));
646  }
647 };
648 
649 } // end namespace internal
650 
654 template <typename MatrixType, typename PermutationIndex>
655 typename ColPivHouseholderQR<MatrixType, PermutationIndex>::HouseholderSequenceType
657  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
658  return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
659 }
660 
665 template <typename Derived>
666 template <typename PermutationIndexType>
670 }
671 
672 } // end namespace Eigen
673 
674 #endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
RealScalar maxPivot() const
Definition: ColPivHouseholderQR.h:398
Index dimensionOfKernel() const
Definition: ColPivHouseholderQR.h:276
ColPivHouseholderQR & setThreshold(Default_t)
Definition: ColPivHouseholderQR.h:366
Index nonzeroPivots() const
Definition: ColPivHouseholderQR.h:390
bool isInjective() const
Definition: ColPivHouseholderQR.h:288
const PermutationType & colsPermutation() const
Definition: ColPivHouseholderQR.h:192
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
Complete orthogonal decomposition (COD) of a matrix.
Definition: ForwardDeclarations.h:433
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
constexpr Index rows() const noexcept
Definition: EigenBase.h:59
MatrixType::RealScalar logAbsDeterminant() const
Definition: ColPivHouseholderQR.h:457
MatrixType::Scalar determinant() const
Definition: ColPivHouseholderQR.h:440
ColPivHouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: ColPivHouseholderQR.h:140
Definition: EigenBase.h:33
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:445
Expression of the inverse of another expression.
Definition: Inverse.h:43
bool isInvertible() const
Definition: ColPivHouseholderQR.h:311
constexpr Index cols() const noexcept
Definition: EigenBase.h:61
const Inverse< ColPivHouseholderQR > inverse() const
Definition: ColPivHouseholderQR.h:321
ComputationInfo info() const
Reports whether the QR factorization was successful.
Definition: ColPivHouseholderQR.h:406
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition: ForwardDeclarations.h:429
const MatrixType & matrixQR() const
Definition: ColPivHouseholderQR.h:169
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Index rank() const
Definition: ColPivHouseholderQR.h:261
ColPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: ColPivHouseholderQR.h:112
bool isSurjective() const
Definition: ColPivHouseholderQR.h:300
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const MatrixType & matrixR() const
Definition: ColPivHouseholderQR.h:183
ColPivHouseholderQR()
Default Constructor.
Definition: ColPivHouseholderQR.h:95
RealScalar threshold() const
Definition: ColPivHouseholderQR.h:375
MatrixType::Scalar signDeterminant() const
Definition: ColPivHouseholderQR.h:464
HouseholderSequenceType householderQ() const
Definition: ColPivHouseholderQR.h:656
const Solve< ColPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Pseudo expression representing a solving operation.
Definition: Solve.h:62
void resize(Index newSize)
Definition: PermutationMatrix.h:122
ComputationInfo
Definition: Constants.h:438
ColPivHouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: ColPivHouseholderQR.h:127
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
MatrixType::RealScalar absDeterminant() const
Definition: ColPivHouseholderQR.h:449
const HCoeffsType & hCoeffs() const
Definition: ColPivHouseholderQR.h:333
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: ColPivHouseholderQR.h:352