Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 <[email protected]>
5 // Copyright (C) 2009 Benoit Jacob <[email protected]>
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 namespace Eigen {
15 
16 namespace internal {
17 
18 template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType;
19 
20 template<typename MatrixType>
21 struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
22 {
23  typedef typename MatrixType::PlainObject ReturnType;
24 };
25 
26 }
27 
49 template<typename _MatrixType> class FullPivHouseholderQR
50 {
51  public:
52 
53  typedef _MatrixType MatrixType;
54  enum {
55  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57  Options = MatrixType::Options,
58  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
59  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
60  };
61  typedef typename MatrixType::Scalar Scalar;
62  typedef typename MatrixType::RealScalar RealScalar;
63  typedef typename MatrixType::Index Index;
64  typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType> MatrixQReturnType;
65  typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
66  typedef Matrix<Index, 1,
67  EIGEN_SIZE_MIN_PREFER_DYNAMIC(ColsAtCompileTime,RowsAtCompileTime), RowMajor, 1,
68  EIGEN_SIZE_MIN_PREFER_FIXED(MaxColsAtCompileTime,MaxRowsAtCompileTime)> IntDiagSizeVectorType;
69  typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
70  typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
71  typedef typename internal::plain_col_type<MatrixType>::type ColVectorType;
72 
79  : m_qr(),
80  m_hCoeffs(),
81  m_rows_transpositions(),
82  m_cols_transpositions(),
83  m_cols_permutation(),
84  m_temp(),
85  m_isInitialized(false),
86  m_usePrescribedThreshold(false) {}
87 
94  FullPivHouseholderQR(Index rows, Index cols)
95  : m_qr(rows, cols),
96  m_hCoeffs((std::min)(rows,cols)),
97  m_rows_transpositions((std::min)(rows,cols)),
98  m_cols_transpositions((std::min)(rows,cols)),
99  m_cols_permutation(cols),
100  m_temp(cols),
101  m_isInitialized(false),
102  m_usePrescribedThreshold(false) {}
103 
116  FullPivHouseholderQR(const MatrixType& matrix)
117  : m_qr(matrix.rows(), matrix.cols()),
118  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
119  m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
120  m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
121  m_cols_permutation(matrix.cols()),
122  m_temp(matrix.cols()),
123  m_isInitialized(false),
124  m_usePrescribedThreshold(false)
125  {
126  compute(matrix);
127  }
128 
147  template<typename Rhs>
148  inline const internal::solve_retval<FullPivHouseholderQR, Rhs>
149  solve(const MatrixBase<Rhs>& b) const
150  {
151  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
152  return internal::solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived());
153  }
154 
157  MatrixQReturnType matrixQ(void) const;
158 
161  const MatrixType& matrixQR() const
162  {
163  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
164  return m_qr;
165  }
166 
167  FullPivHouseholderQR& compute(const MatrixType& matrix);
168 
171  {
172  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
173  return m_cols_permutation;
174  }
175 
178  {
179  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
180  return m_rows_transpositions;
181  }
182 
196  typename MatrixType::RealScalar absDeterminant() const;
197 
210  typename MatrixType::RealScalar logAbsDeterminant() const;
211 
218  inline Index rank() const
219  {
220  using std::abs;
221  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
222  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
223  Index result = 0;
224  for(Index i = 0; i < m_nonzero_pivots; ++i)
225  result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
226  return result;
227  }
228 
235  inline Index dimensionOfKernel() const
236  {
237  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
238  return cols() - rank();
239  }
240 
248  inline bool isInjective() const
249  {
250  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
251  return rank() == cols();
252  }
253 
261  inline bool isSurjective() const
262  {
263  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
264  return rank() == rows();
265  }
266 
273  inline bool isInvertible() const
274  {
275  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
276  return isInjective() && isSurjective();
277  }
278  inline const
284  internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType>
285  inverse() const
286  {
287  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
288  return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType>
289  (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
290  }
291 
292  inline Index rows() const { return m_qr.rows(); }
293  inline Index cols() const { return m_qr.cols(); }
294 
299  const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
300 
319  {
320  m_usePrescribedThreshold = true;
321  m_prescribedThreshold = threshold;
322  return *this;
323  }
324 
334  {
335  m_usePrescribedThreshold = false;
336  return *this;
337  }
338 
343  RealScalar threshold() const
344  {
345  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
346  return m_usePrescribedThreshold ? m_prescribedThreshold
347  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
348  // and turns out to be identical to Higham's formula used already in LDLt.
349  : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize());
350  }
351 
359  inline Index nonzeroPivots() const
360  {
361  eigen_assert(m_isInitialized && "LU is not initialized.");
362  return m_nonzero_pivots;
363  }
364 
368  RealScalar maxPivot() const { return m_maxpivot; }
369 
370  protected:
371 
372  static void check_template_parameters()
373  {
374  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
375  }
376 
377  MatrixType m_qr;
378  HCoeffsType m_hCoeffs;
379  IntDiagSizeVectorType m_rows_transpositions;
380  IntDiagSizeVectorType m_cols_transpositions;
381  PermutationType m_cols_permutation;
382  RowVectorType m_temp;
383  bool m_isInitialized, m_usePrescribedThreshold;
384  RealScalar m_prescribedThreshold, m_maxpivot;
385  Index m_nonzero_pivots;
386  RealScalar m_precision;
387  Index m_det_pq;
388 };
389 
390 template<typename MatrixType>
391 typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::absDeterminant() const
392 {
393  using std::abs;
394  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
395  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
396  return abs(m_qr.diagonal().prod());
397 }
398 
399 template<typename MatrixType>
400 typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::logAbsDeterminant() const
401 {
402  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
403  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
404  return m_qr.diagonal().cwiseAbs().array().log().sum();
405 }
406 
413 template<typename MatrixType>
415 {
416  check_template_parameters();
417 
418  using std::abs;
419  Index rows = matrix.rows();
420  Index cols = matrix.cols();
421  Index size = (std::min)(rows,cols);
422 
423  m_qr = matrix;
424  m_hCoeffs.resize(size);
425 
426  m_temp.resize(cols);
427 
428  m_precision = NumTraits<Scalar>::epsilon() * RealScalar(size);
429 
430  m_rows_transpositions.resize(size);
431  m_cols_transpositions.resize(size);
432  Index number_of_transpositions = 0;
433 
434  RealScalar biggest(0);
435 
436  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
437  m_maxpivot = RealScalar(0);
438 
439  for (Index k = 0; k < size; ++k)
440  {
441  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
442  RealScalar biggest_in_corner;
443 
444  biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
445  .cwiseAbs()
446  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
447  row_of_biggest_in_corner += k;
448  col_of_biggest_in_corner += k;
449  if(k==0) biggest = biggest_in_corner;
450 
451  // if the corner is negligible, then we have less than full rank, and we can finish early
452  if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
453  {
454  m_nonzero_pivots = k;
455  for(Index i = k; i < size; i++)
456  {
457  m_rows_transpositions.coeffRef(i) = i;
458  m_cols_transpositions.coeffRef(i) = i;
459  m_hCoeffs.coeffRef(i) = Scalar(0);
460  }
461  break;
462  }
463 
464  m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
465  m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
466  if(k != row_of_biggest_in_corner) {
467  m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
468  ++number_of_transpositions;
469  }
470  if(k != col_of_biggest_in_corner) {
471  m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
472  ++number_of_transpositions;
473  }
474 
475  RealScalar beta;
476  m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
477  m_qr.coeffRef(k,k) = beta;
478 
479  // remember the maximum absolute value of diagonal coefficients
480  if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
481 
482  m_qr.bottomRightCorner(rows-k, cols-k-1)
483  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
484  }
485 
486  m_cols_permutation.setIdentity(cols);
487  for(Index k = 0; k < size; ++k)
488  m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
489 
490  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
491  m_isInitialized = true;
492 
493  return *this;
494 }
495 
496 namespace internal {
497 
498 template<typename _MatrixType, typename Rhs>
499 struct solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
500  : solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs>
501 {
502  EIGEN_MAKE_SOLVE_HELPERS(FullPivHouseholderQR<_MatrixType>,Rhs)
503 
504  template<typename Dest> void evalTo(Dest& dst) const
505  {
506  const Index rows = dec().rows(), cols = dec().cols();
507  eigen_assert(rhs().rows() == rows);
508 
509  // FIXME introduce nonzeroPivots() and use it here. and more generally,
510  // make the same improvements in this dec as in FullPivLU.
511  if(dec().rank()==0)
512  {
513  dst.setZero();
514  return;
515  }
516 
517  typename Rhs::PlainObject c(rhs());
518 
519  Matrix<Scalar,1,Rhs::ColsAtCompileTime> temp(rhs().cols());
520  for (Index k = 0; k < dec().rank(); ++k)
521  {
522  Index remainingSize = rows-k;
523  c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k)));
524  c.bottomRightCorner(remainingSize, rhs().cols())
525  .applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1),
526  dec().hCoeffs().coeff(k), &temp.coeffRef(0));
527  }
528 
529  dec().matrixQR()
530  .topLeftCorner(dec().rank(), dec().rank())
531  .template triangularView<Upper>()
532  .solveInPlace(c.topRows(dec().rank()));
533 
534  for(Index i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
535  for(Index i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
536  }
537 };
538 
545 template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType
546  : public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
547 {
548 public:
549  typedef typename MatrixType::Index Index;
550  typedef typename FullPivHouseholderQR<MatrixType>::IntDiagSizeVectorType IntDiagSizeVectorType;
551  typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
552  typedef Matrix<typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime, RowMajor, 1,
553  MatrixType::MaxRowsAtCompileTime> WorkVectorType;
554 
555  FullPivHouseholderQRMatrixQReturnType(const MatrixType& qr,
556  const HCoeffsType& hCoeffs,
557  const IntDiagSizeVectorType& rowsTranspositions)
558  : m_qr(qr),
559  m_hCoeffs(hCoeffs),
560  m_rowsTranspositions(rowsTranspositions)
561  {}
562 
563  template <typename ResultType>
564  void evalTo(ResultType& result) const
565  {
566  const Index rows = m_qr.rows();
567  WorkVectorType workspace(rows);
568  evalTo(result, workspace);
569  }
570 
571  template <typename ResultType>
572  void evalTo(ResultType& result, WorkVectorType& workspace) const
573  {
574  using numext::conj;
575  // compute the product H'_0 H'_1 ... H'_n-1,
576  // where H_k is the k-th Householder transformation I - h_k v_k v_k'
577  // and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
578  const Index rows = m_qr.rows();
579  const Index cols = m_qr.cols();
580  const Index size = (std::min)(rows, cols);
581  workspace.resize(rows);
582  result.setIdentity(rows, rows);
583  for (Index k = size-1; k >= 0; k--)
584  {
585  result.block(k, k, rows-k, rows-k)
586  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k));
587  result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
588  }
589  }
590 
591  Index rows() const { return m_qr.rows(); }
592  Index cols() const { return m_qr.rows(); }
593 
594 protected:
595  typename MatrixType::Nested m_qr;
596  typename HCoeffsType::Nested m_hCoeffs;
597  typename IntDiagSizeVectorType::Nested m_rowsTranspositions;
598 };
599 
600 } // end namespace internal
601 
602 template<typename MatrixType>
603 inline typename FullPivHouseholderQR<MatrixType>::MatrixQReturnType FullPivHouseholderQR<MatrixType>::matrixQ() const
604 {
605  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
606  return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
607 }
608 
613 template<typename Derived>
616 {
617  return FullPivHouseholderQR<PlainObject>(eval());
618 }
619 
620 } // end namespace Eigen
621 
622 #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
MatrixType::RealScalar absDeterminant() const
Definition: FullPivHouseholderQR.h:391
const internal::solve_retval< FullPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivHouseholderQR.h:149
bool isInjective() const
Definition: FullPivHouseholderQR.h:248
Householder rank-revealing QR decomposition of a matrix with full pivoting.
Definition: ForwardDeclarations.h:223
const PermutationType & colsPermutation() const
Definition: FullPivHouseholderQR.h:170
Index dimensionOfKernel() const
Definition: FullPivHouseholderQR.h:235
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const HCoeffsType & hCoeffs() const
Definition: FullPivHouseholderQR.h:299
FullPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: FullPivHouseholderQR.h:94
bool isInvertible() const
Definition: FullPivHouseholderQR.h:273
bool isSurjective() const
Definition: FullPivHouseholderQR.h:261
FullPivHouseholderQR & compute(const MatrixType &matrix)
Definition: FullPivHouseholderQR.h:414
RealScalar maxPivot() const
Definition: FullPivHouseholderQR.h:368
const FullPivHouseholderQR< PlainObject > fullPivHouseholderQr() const
Definition: FullPivHouseholderQR.h:615
FullPivHouseholderQR()
Default Constructor.
Definition: FullPivHouseholderQR.h:78
Index rank() const
Definition: FullPivHouseholderQR.h:218
RealScalar threshold() const
Definition: FullPivHouseholderQR.h:343
const MatrixType & matrixQR() const
Definition: FullPivHouseholderQR.h:161
FullPivHouseholderQR & setThreshold(Default_t)
Definition: FullPivHouseholderQR.h:333
FullPivHouseholderQR(const MatrixType &matrix)
Constructs a QR factorization from a given matrix.
Definition: FullPivHouseholderQR.h:116
Index nonzeroPivots() const
Definition: FullPivHouseholderQR.h:359
const internal::solve_retval< FullPivHouseholderQR, typename MatrixType::IdentityReturnType > inverse() const
Definition: FullPivHouseholderQR.h:285
MatrixType::RealScalar logAbsDeterminant() const
Definition: FullPivHouseholderQR.h:400
Definition: Constants.h:266
const IntDiagSizeVectorType & rowsTranspositions() const
Definition: FullPivHouseholderQR.h:177
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
MatrixQReturnType matrixQ(void) const
Definition: FullPivHouseholderQR.h:603
FullPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: FullPivHouseholderQR.h:318