$darkmode
Eigen  5.0.1-dev
HouseholderQR.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
6 // Copyright (C) 2010 Vincent Lejeune
7 //
8 // This Source Code Form is subject to the terms of the Mozilla
9 // Public License v. 2.0. If a copy of the MPL was not distributed
10 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11 
12 #ifndef EIGEN_QR_H
13 #define EIGEN_QR_H
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21 template <typename MatrixType_>
22 struct traits<HouseholderQR<MatrixType_>> : traits<MatrixType_> {
23  typedef MatrixXpr XprKind;
24  typedef SolverStorage StorageKind;
25  typedef int StorageIndex;
26  enum { Flags = 0 };
27 };
28 
29 } // end namespace internal
30 
58 template <typename MatrixType_>
59 class HouseholderQR : public SolverBase<HouseholderQR<MatrixType_>> {
60  public:
61  typedef MatrixType_ MatrixType;
62  typedef SolverBase<HouseholderQR> Base;
63  friend class SolverBase<HouseholderQR>;
64 
65  EIGEN_GENERIC_PUBLIC_INTERFACE(HouseholderQR)
66  enum {
67  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
68  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
69  };
70  typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, (MatrixType::Flags & RowMajorBit) ? RowMajor : ColMajor,
71  MaxRowsAtCompileTime, MaxRowsAtCompileTime>
72  MatrixQType;
73  typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
74  typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
75  typedef HouseholderSequence<MatrixType, internal::remove_all_t<typename HCoeffsType::ConjugateReturnType>>
76  HouseholderSequenceType;
77 
85  eigen_assert(m_isInitialized && "HouseHolderQR is not initialized.");
86  return Success;
87  }
88 
95  HouseholderQR() : m_qr(), m_hCoeffs(), m_temp(), m_isInitialized(false) {}
96 
104  : m_qr(rows, cols), m_hCoeffs((std::min)(rows, cols)), m_temp(cols), m_isInitialized(false) {}
105 
118  template <typename InputType>
119  explicit HouseholderQR(const EigenBase<InputType>& matrix)
120  : m_qr(matrix.rows(), matrix.cols()),
121  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
122  m_temp(matrix.cols()),
123  m_isInitialized(false) {
124  compute(matrix.derived());
125  }
126 
134  template <typename InputType>
136  : m_qr(matrix.derived()),
137  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
138  m_temp(matrix.cols()),
139  m_isInitialized(false) {
140  computeInPlace();
141  }
142 
143 #ifdef EIGEN_PARSED_BY_DOXYGEN
144 
158  template <typename Rhs>
159  inline const Solve<HouseholderQR, Rhs> solve(const MatrixBase<Rhs>& b) const;
160 #endif
161 
172  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
173  return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
174  }
175 
179  const MatrixType& matrixQR() const {
180  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
181  return m_qr;
182  }
183 
184  template <typename InputType>
185  HouseholderQR& compute(const EigenBase<InputType>& matrix) {
186  m_qr = matrix.derived();
187  computeInPlace();
188  return *this;
189  }
190 
206  typename MatrixType::Scalar determinant() const;
207 
223  typename MatrixType::RealScalar absDeterminant() const;
224 
240  typename MatrixType::RealScalar logAbsDeterminant() const;
241 
257  typename MatrixType::Scalar signDeterminant() const;
258 
259  inline Index rows() const { return m_qr.rows(); }
260  inline Index cols() const { return m_qr.cols(); }
261 
266  const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
267 
268 #ifndef EIGEN_PARSED_BY_DOXYGEN
269  template <typename RhsType, typename DstType>
270  void _solve_impl(const RhsType& rhs, DstType& dst) const;
271 
272  template <bool Conjugate, typename RhsType, typename DstType>
273  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
274 #endif
275 
276  protected:
277  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
278 
279  void computeInPlace();
280 
281  MatrixType m_qr;
282  HCoeffsType m_hCoeffs;
283  RowVectorType m_temp;
284  bool m_isInitialized;
285 };
286 
287 namespace internal {
288 
290 template <typename HCoeffs, typename Scalar, bool IsComplex>
291 struct householder_determinant {
292  static void run(const HCoeffs& hCoeffs, Scalar& out_det) {
293  out_det = Scalar(1);
294  Index size = hCoeffs.rows();
295  for (Index i = 0; i < size; i++) {
296  // For each valid reflection Q_n,
297  // det(Q_n) = - conj(h_n) / h_n
298  // where h_n is the Householder coefficient.
299  if (hCoeffs(i) != Scalar(0)) out_det *= -numext::conj(hCoeffs(i)) / hCoeffs(i);
300  }
301  }
302 };
303 
305 template <typename HCoeffs, typename Scalar>
306 struct householder_determinant<HCoeffs, Scalar, false> {
307  static void run(const HCoeffs& hCoeffs, Scalar& out_det) {
308  bool negated = false;
309  Index size = hCoeffs.rows();
310  for (Index i = 0; i < size; i++) {
311  // Each valid reflection negates the determinant.
312  if (hCoeffs(i) != Scalar(0)) negated ^= true;
313  }
314  out_det = negated ? Scalar(-1) : Scalar(1);
315  }
316 };
317 
318 } // end namespace internal
319 
320 template <typename MatrixType>
321 typename MatrixType::Scalar HouseholderQR<MatrixType>::determinant() const {
322  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
323  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
324  Scalar detQ;
325  internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
326  return m_qr.diagonal().prod() * detQ;
327 }
328 
329 template <typename MatrixType>
330 typename MatrixType::RealScalar HouseholderQR<MatrixType>::absDeterminant() const {
331  using std::abs;
332  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
333  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
334  return abs(m_qr.diagonal().prod());
335 }
336 
337 template <typename MatrixType>
338 typename MatrixType::RealScalar HouseholderQR<MatrixType>::logAbsDeterminant() const {
339  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
340  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
341  return m_qr.diagonal().cwiseAbs().array().log().sum();
342 }
343 
344 template <typename MatrixType>
345 typename MatrixType::Scalar HouseholderQR<MatrixType>::signDeterminant() const {
346  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
347  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
348  Scalar detQ;
349  internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
350  return detQ * m_qr.diagonal().array().sign().prod();
351 }
352 
353 namespace internal {
354 
356 template <typename MatrixQR, typename HCoeffs>
357 void householder_qr_inplace_unblocked(MatrixQR& mat, HCoeffs& hCoeffs, typename MatrixQR::Scalar* tempData = 0) {
358  typedef typename MatrixQR::Scalar Scalar;
359  typedef typename MatrixQR::RealScalar RealScalar;
360  Index rows = mat.rows();
361  Index cols = mat.cols();
362  Index size = (std::min)(rows, cols);
363 
364  eigen_assert(hCoeffs.size() == size);
365 
367  TempType tempVector;
368  if (tempData == 0) {
369  tempVector.resize(cols);
370  tempData = tempVector.data();
371  }
372 
373  for (Index k = 0; k < size; ++k) {
374  Index remainingRows = rows - k;
375  Index remainingCols = cols - k - 1;
376 
377  RealScalar beta;
378  mat.col(k).tail(remainingRows).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
379  mat.coeffRef(k, k) = beta;
380 
381  // apply H to remaining part of m_qr from the left
382  mat.bottomRightCorner(remainingRows, remainingCols)
383  .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows - 1), hCoeffs.coeffRef(k), tempData + k + 1);
384  }
385 }
386 
387 // TODO: add a corresponding public API for updating a QR factorization
396 template <typename MatrixQR, typename HCoeffs, typename VectorQR>
397 void householder_qr_inplace_update(MatrixQR& mat, HCoeffs& hCoeffs, const VectorQR& newColumn,
398  typename MatrixQR::Index k, typename MatrixQR::Scalar* tempData) {
399  typedef typename MatrixQR::Index Index;
400  typedef typename MatrixQR::RealScalar RealScalar;
401  Index rows = mat.rows();
402 
403  eigen_assert(k < mat.cols());
404  eigen_assert(k < rows);
405  eigen_assert(hCoeffs.size() == mat.cols());
406  eigen_assert(newColumn.size() == rows);
407  eigen_assert(tempData);
408 
409  // Store new column in mat at column k
410  mat.col(k) = newColumn;
411  // Apply H = H_1...H_{k-1} on newColumn (skip if k=0)
412  for (Index i = 0; i < k; ++i) {
413  Index remainingRows = rows - i;
414  mat.col(k)
415  .tail(remainingRows)
416  .applyHouseholderOnTheLeft(mat.col(i).tail(remainingRows - 1), hCoeffs.coeffRef(i), tempData + i + 1);
417  }
418  // Construct Householder projector in-place in column k
419  RealScalar beta;
420  mat.col(k).tail(rows - k).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
421  mat.coeffRef(k, k) = beta;
422 }
423 
425 template <typename MatrixQR, typename HCoeffs, typename MatrixQRScalar = typename MatrixQR::Scalar,
426  bool InnerStrideIsOne = (MatrixQR::InnerStrideAtCompileTime == 1 && HCoeffs::InnerStrideAtCompileTime == 1)>
427 struct householder_qr_inplace_blocked {
428  // This is specialized for LAPACK-supported Scalar types in HouseholderQR_LAPACKE.h
429  static void run(MatrixQR& mat, HCoeffs& hCoeffs, Index maxBlockSize = 32, typename MatrixQR::Scalar* tempData = 0) {
430  typedef typename MatrixQR::Scalar Scalar;
431  typedef Block<MatrixQR, Dynamic, Dynamic> BlockType;
432 
433  Index rows = mat.rows();
434  Index cols = mat.cols();
435  Index size = (std::min)(rows, cols);
436 
437  typedef Matrix<Scalar, Dynamic, 1, ColMajor, MatrixQR::MaxColsAtCompileTime, 1> TempType;
438  TempType tempVector;
439  if (tempData == 0) {
440  tempVector.resize(cols);
441  tempData = tempVector.data();
442  }
443 
444  Index blockSize = (std::min)(maxBlockSize, size);
445 
446  Index k = 0;
447  for (k = 0; k < size; k += blockSize) {
448  Index bs = (std::min)(size - k, blockSize); // actual size of the block
449  Index tcols = cols - k - bs; // trailing columns
450  Index brows = rows - k; // rows of the block
451 
452  // partition the matrix:
453  // A00 | A01 | A02
454  // mat = A10 | A11 | A12
455  // A20 | A21 | A22
456  // and performs the qr dec of [A11^T A12^T]^T
457  // and update [A21^T A22^T]^T using level 3 operations.
458  // Finally, the algorithm continue on A22
459 
460  BlockType A11_21 = mat.block(k, k, brows, bs);
461  Block<HCoeffs, Dynamic, 1> hCoeffsSegment = hCoeffs.segment(k, bs);
462 
463  householder_qr_inplace_unblocked(A11_21, hCoeffsSegment, tempData);
464 
465  if (tcols) {
466  BlockType A21_22 = mat.block(k, k + bs, brows, tcols);
467  apply_block_householder_on_the_left(A21_22, A11_21, hCoeffsSegment, false); // false == backward
468  }
469  }
470  }
471 };
472 
473 } // end namespace internal
474 
475 #ifndef EIGEN_PARSED_BY_DOXYGEN
476 template <typename MatrixType_>
477 template <typename RhsType, typename DstType>
478 void HouseholderQR<MatrixType_>::_solve_impl(const RhsType& rhs, DstType& dst) const {
479  const Index rank = (std::min)(rows(), cols());
480 
481  typename RhsType::PlainObject c(rhs);
482 
483  c.applyOnTheLeft(householderQ().setLength(rank).adjoint());
484 
485  m_qr.topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(c.topRows(rank));
486 
487  dst.topRows(rank) = c.topRows(rank);
488  dst.bottomRows(cols() - rank).setZero();
489 }
490 
491 template <typename MatrixType_>
492 template <bool Conjugate, typename RhsType, typename DstType>
493 void HouseholderQR<MatrixType_>::_solve_impl_transposed(const RhsType& rhs, DstType& dst) const {
494  const Index rank = (std::min)(rows(), cols());
495 
496  typename RhsType::PlainObject c(rhs);
497 
498  m_qr.topLeftCorner(rank, rank)
499  .template triangularView<Upper>()
500  .transpose()
501  .template conjugateIf<Conjugate>()
502  .solveInPlace(c.topRows(rank));
503 
504  dst.topRows(rank) = c.topRows(rank);
505  dst.bottomRows(rows() - rank).setZero();
506 
507  dst.applyOnTheLeft(householderQ().setLength(rank).template conjugateIf<!Conjugate>());
508 }
509 #endif
510 
517 template <typename MatrixType>
519  Index rows = m_qr.rows();
520  Index cols = m_qr.cols();
521  Index size = (std::min)(rows, cols);
522 
523  m_hCoeffs.resize(size);
524 
525  m_temp.resize(cols);
526 
527  internal::householder_qr_inplace_blocked<MatrixType, HCoeffsType>::run(m_qr, m_hCoeffs, 48, m_temp.data());
528 
529  m_isInitialized = true;
530 }
531 
536 template <typename Derived>
538  return HouseholderQR<PlainObject>(eval());
539 }
540 
541 } // end namespace Eigen
542 
543 #endif // EIGEN_QR_H
HouseholderSequenceType householderQ() const
Definition: HouseholderQR.h:171
const MatrixType & matrixQR() const
Definition: HouseholderQR.h:179
Definition: Constants.h:318
constexpr void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:268
MatrixType::Scalar signDeterminant() const
Definition: HouseholderQR.h:345
MatrixType::RealScalar absDeterminant() const
Definition: HouseholderQR.h:330
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
MatrixType::Scalar determinant() const
Definition: HouseholderQR.h:321
Definition: BFloat16.h:231
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
const unsigned int RowMajorBit
Definition: Constants.h:70
Definition: EigenBase.h:33
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:445
HouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: HouseholderQR.h:135
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
MatrixType::RealScalar logAbsDeterminant() const
Definition: HouseholderQR.h:338
const Solve< HouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
ComputationInfo info() const
Reports whether the QR factorization was successful.
Definition: HouseholderQR.h:84
const HCoeffsType & hCoeffs() const
Definition: HouseholderQR.h:266
Householder QR decomposition of a matrix.
Definition: ForwardDeclarations.h:427
HouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: HouseholderQR.h:103
void computeInPlace()
Definition: HouseholderQR.h:518
Definition: Constants.h:320
Pseudo expression representing a solving operation.
Definition: Solve.h:62
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:186
ComputationInfo
Definition: Constants.h:438
HouseholderQR()
Default Constructor.
Definition: HouseholderQR.h:95
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
HouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: HouseholderQR.h:119