$darkmode
Eigen  5.0.1-dev
FullPivHouseholderQR.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_FULLPIVOTINGHOUSEHOLDERQR_H
12 #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 
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;
26  enum { Flags = 0 };
27 };
28 
29 template <typename MatrixType, typename PermutationIndex>
31 
32 template <typename MatrixType, typename PermutationIndex>
33 struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> > {
34  typedef typename MatrixType::PlainObject ReturnType;
35 };
36 
37 } // end namespace internal
38 
62 template <typename MatrixType_, typename PermutationIndex_>
63 class FullPivHouseholderQR : public SolverBase<FullPivHouseholderQR<MatrixType_, PermutationIndex_> > {
64  public:
65  typedef MatrixType_ MatrixType;
67  friend class SolverBase<FullPivHouseholderQR>;
68  typedef PermutationIndex_ PermutationIndex;
69  EIGEN_GENERIC_PUBLIC_INTERFACE(FullPivHouseholderQR)
70 
71  enum {
72  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
74  };
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;
84 
92  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
93  return Success;
94  }
95 
102  : m_qr(),
103  m_hCoeffs(),
104  m_rows_transpositions(),
105  m_cols_transpositions(),
106  m_cols_permutation(),
107  m_temp(),
108  m_isInitialized(false),
109  m_usePrescribedThreshold(false) {}
110 
118  : m_qr(rows, cols),
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),
123  m_temp(cols),
124  m_isInitialized(false),
125  m_usePrescribedThreshold(false) {}
126 
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) {
149  compute(matrix.derived());
150  }
151 
159  template <typename InputType>
161  : m_qr(matrix.derived()),
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) {
169  computeInPlace();
170  }
171 
172 #ifdef EIGEN_PARSED_BY_DOXYGEN
173 
188  template <typename Rhs>
189  inline const Solve<FullPivHouseholderQR, Rhs> solve(const MatrixBase<Rhs>& b) const;
190 #endif
191 
194  MatrixQReturnType matrixQ(void) const;
195 
198  const MatrixType& matrixQR() const {
199  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
200  return m_qr;
201  }
202 
203  template <typename InputType>
204  FullPivHouseholderQR& compute(const EigenBase<InputType>& matrix);
205 
208  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
209  return m_cols_permutation;
210  }
211 
214  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
215  return m_rows_transpositions;
216  }
217 
231  typename MatrixType::Scalar determinant() const;
232 
246  typename MatrixType::RealScalar absDeterminant() const;
247 
260  typename MatrixType::RealScalar logAbsDeterminant() const;
261 
274  typename MatrixType::Scalar signDeterminant() const;
275 
282  inline Index rank() const {
283  using std::abs;
284  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
285  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
286  Index result = 0;
287  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_qr.coeff(i, i)) > premultiplied_threshold);
288  return result;
289  }
290 
297  inline Index dimensionOfKernel() const {
298  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
299  return cols() - rank();
300  }
301 
309  inline bool isInjective() const {
310  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
311  return rank() == cols();
312  }
313 
321  inline bool isSurjective() const {
322  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
323  return rank() == rows();
324  }
325 
332  inline bool isInvertible() const {
333  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
334  return isInjective() && isSurjective();
335  }
336 
342  inline const Inverse<FullPivHouseholderQR> inverse() const {
343  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
344  return Inverse<FullPivHouseholderQR>(*this);
345  }
346 
347  inline Index rows() const { return m_qr.rows(); }
348  inline Index cols() const { return m_qr.cols(); }
349 
354  const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
355 
374  m_usePrescribedThreshold = true;
375  m_prescribedThreshold = threshold;
376  return *this;
377  }
378 
388  m_usePrescribedThreshold = false;
389  return *this;
390  }
391 
396  RealScalar threshold() const {
397  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
398  return m_usePrescribedThreshold ? m_prescribedThreshold
399  // this formula comes from experimenting (see "LU precision tuning" thread on the
400  // list) and turns out to be identical to Higham's formula used already in LDLt.
401  : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize());
402  }
403 
411  inline Index nonzeroPivots() const {
412  eigen_assert(m_isInitialized && "LU is not initialized.");
413  return m_nonzero_pivots;
414  }
415 
419  RealScalar maxPivot() const { return m_maxpivot; }
420 
421 #ifndef EIGEN_PARSED_BY_DOXYGEN
422  template <typename RhsType, typename DstType>
423  void _solve_impl(const RhsType& rhs, DstType& dst) const;
424 
425  template <bool Conjugate, typename RhsType, typename DstType>
426  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
427 #endif
428 
429  protected:
430  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
431 
432  void computeInPlace();
433 
434  MatrixType m_qr;
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;
444  Index m_det_p;
445 };
446 
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!");
451  Scalar detQ;
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);
454 }
455 
456 template <typename MatrixType, typename PermutationIndex>
458  using std::abs;
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);
462 }
463 
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!");
468  return isInjective() ? m_qr.diagonal().cwiseAbs().array().log().sum() : -NumTraits<RealScalar>::infinity();
469 }
470 
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!");
475  Scalar detQ;
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);
478 }
479 
486 template <typename MatrixType, typename PermutationIndex>
487 template <typename InputType>
489  const EigenBase<InputType>& matrix) {
490  m_qr = matrix.derived();
491  computeInPlace();
492  return *this;
493 }
494 
495 template <typename MatrixType, typename PermutationIndex>
497  eigen_assert(m_qr.cols() <= NumTraits<PermutationIndex>::highest());
498  using std::abs;
499  Index rows = m_qr.rows();
500  Index cols = m_qr.cols();
501  Index size = (std::min)(rows, cols);
502 
503  m_hCoeffs.resize(size);
504 
505  m_temp.resize(cols);
506 
507  m_precision = NumTraits<Scalar>::epsilon() * RealScalar(size);
508 
509  m_rows_transpositions.resize(size);
510  m_cols_transpositions.resize(size);
511  Index number_of_transpositions = 0;
512 
513  RealScalar biggest(0);
514 
515  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
516  m_maxpivot = RealScalar(0);
517 
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;
522 
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;
531 
532  // if the corner is negligible, then we have less than full rank, and we can finish early
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);
539  }
540  break;
541  }
542 
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;
548  }
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;
552  }
553 
554  RealScalar beta;
555  m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
556  m_qr.coeffRef(k, k) = beta;
557 
558  // remember the maximum absolute value of diagonal coefficients
559  if (abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
560 
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));
563  }
564 
565  m_cols_permutation.setIdentity(cols);
566  for (Index k = 0; k < size; ++k) m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
567 
568  m_det_p = (number_of_transpositions % 2) ? -1 : 1;
569  m_isInitialized = true;
570 }
571 
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();
577 
578  // FIXME introduce nonzeroPivots() and use it here. and more generally,
579  // make the same improvements in this dec as in FullPivLU.
580  if (l_rank == 0) {
581  dst.setZero();
582  return;
583  }
584 
585  typename RhsType::PlainObject c(rhs);
586 
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));
593  }
594 
595  m_qr.topLeftCorner(l_rank, l_rank).template triangularView<Upper>().solveInPlace(c.topRows(l_rank));
596 
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();
599 }
600 
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();
606 
607  if (l_rank == 0) {
608  dst.setZero();
609  return;
610  }
611 
612  typename RhsType::PlainObject c(m_cols_permutation.transpose() * rhs);
613 
614  m_qr.topLeftCorner(l_rank, l_rank)
615  .template triangularView<Upper>()
616  .transpose()
617  .template conjugateIf<Conjugate>()
618  .solveInPlace(c.topRows(l_rank));
619 
620  dst.topRows(l_rank) = c.topRows(l_rank);
621  dst.bottomRows(rows() - l_rank).setZero();
622 
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;
627 
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));
631 
632  dst.row(k).swap(dst.row(m_rows_transpositions.coeff(k)));
633  }
634 }
635 #endif
636 
637 namespace internal {
638 
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>,
643  Dense2Dense> {
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()));
649  }
650 };
651 
658 template <typename MatrixType, typename PermutationIndex>
659 struct FullPivHouseholderQRMatrixQReturnType
660  : public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> > {
661  public:
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>
666  WorkVectorType;
667 
668  FullPivHouseholderQRMatrixQReturnType(const MatrixType& qr, const HCoeffsType& hCoeffs,
669  const IntDiagSizeVectorType& rowsTranspositions)
670  : m_qr(qr), m_hCoeffs(hCoeffs), m_rowsTranspositions(rowsTranspositions) {}
671 
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);
677  }
678 
679  template <typename ResultType>
680  void evalTo(ResultType& result, WorkVectorType& workspace) const {
681  using numext::conj;
682  // compute the product H'_0 H'_1 ... H'_n-1,
683  // where H_k is the k-th Householder transformation I - h_k v_k v_k'
684  // and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
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)));
694  }
695  }
696 
697  Index rows() const { return m_qr.rows(); }
698  Index cols() const { return m_qr.rows(); }
699 
700  protected:
701  typename MatrixType::Nested m_qr;
702  typename HCoeffsType::Nested m_hCoeffs;
703  typename IntDiagSizeVectorType::Nested m_rowsTranspositions;
704 };
705 
706 // template<typename MatrixType>
707 // struct evaluator<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
708 // : public evaluator<ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> > >
709 // {};
710 
711 } // end namespace internal
712 
713 template <typename MatrixType, typename PermutationIndex>
714 inline typename FullPivHouseholderQR<MatrixType, PermutationIndex>::MatrixQReturnType
716  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
717  return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
718 }
719 
724 template <typename Derived>
725 template <typename PermutationIndex>
729 }
730 
731 } // end namespace Eigen
732 
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