ColPivHouseholderQR.h
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
00006 //
00007 // This Source Code Form is subject to the terms of the Mozilla
00008 // Public License v. 2.0. If a copy of the MPL was not distributed
00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00010 
00011 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
00012 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
00013 
00014 namespace Eigen { 
00015 
00037 template<typename _MatrixType> class ColPivHouseholderQR
00038 {
00039   public:
00040 
00041     typedef _MatrixType MatrixType;
00042     enum {
00043       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00044       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00045       Options = MatrixType::Options,
00046       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00047       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00048     };
00049     typedef typename MatrixType::Scalar Scalar;
00050     typedef typename MatrixType::RealScalar RealScalar;
00051     typedef typename MatrixType::Index Index;
00052     typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
00053     typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
00054     typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
00055     typedef typename internal::plain_row_type<MatrixType, Index>::type IntRowVectorType;
00056     typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
00057     typedef typename internal::plain_row_type<MatrixType, RealScalar>::type RealRowVectorType;
00058     typedef typename HouseholderSequence<MatrixType,HCoeffsType>::ConjugateReturnType HouseholderSequenceType;
00059     
00060   private:
00061     
00062     typedef typename PermutationType::Index PermIndexType;
00063     
00064   public:
00065 
00072     ColPivHouseholderQR()
00073       : m_qr(),
00074         m_hCoeffs(),
00075         m_colsPermutation(),
00076         m_colsTranspositions(),
00077         m_temp(),
00078         m_colSqNorms(),
00079         m_isInitialized(false) {}
00080 
00087     ColPivHouseholderQR(Index rows, Index cols)
00088       : m_qr(rows, cols),
00089         m_hCoeffs((std::min)(rows,cols)),
00090         m_colsPermutation(PermIndexType(cols)),
00091         m_colsTranspositions(cols),
00092         m_temp(cols),
00093         m_colSqNorms(cols),
00094         m_isInitialized(false),
00095         m_usePrescribedThreshold(false) {}
00096 
00097     ColPivHouseholderQR(const MatrixType& matrix)
00098       : m_qr(matrix.rows(), matrix.cols()),
00099         m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
00100         m_colsPermutation(PermIndexType(matrix.cols())),
00101         m_colsTranspositions(matrix.cols()),
00102         m_temp(matrix.cols()),
00103         m_colSqNorms(matrix.cols()),
00104         m_isInitialized(false),
00105         m_usePrescribedThreshold(false)
00106     {
00107       compute(matrix);
00108     }
00109 
00127     template<typename Rhs>
00128     inline const internal::solve_retval<ColPivHouseholderQR, Rhs>
00129     solve(const MatrixBase<Rhs>& b) const
00130     {
00131       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00132       return internal::solve_retval<ColPivHouseholderQR, Rhs>(*this, b.derived());
00133     }
00134 
00135     HouseholderSequenceType householderQ(void) const;
00136 
00139     const MatrixType& matrixQR() const
00140     {
00141       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00142       return m_qr;
00143     }
00144 
00145     ColPivHouseholderQR& compute(const MatrixType& matrix);
00146 
00147     const PermutationType& colsPermutation() const
00148     {
00149       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00150       return m_colsPermutation;
00151     }
00152 
00166     typename MatrixType::RealScalar absDeterminant() const;
00167 
00180     typename MatrixType::RealScalar logAbsDeterminant() const;
00181 
00188     inline Index rank() const
00189     {
00190       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00191       RealScalar premultiplied_threshold = internal::abs(m_maxpivot) * threshold();
00192       Index result = 0;
00193       for(Index i = 0; i < m_nonzero_pivots; ++i)
00194         result += (internal::abs(m_qr.coeff(i,i)) > premultiplied_threshold);
00195       return result;
00196     }
00197 
00204     inline Index dimensionOfKernel() const
00205     {
00206       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00207       return cols() - rank();
00208     }
00209 
00217     inline bool isInjective() const
00218     {
00219       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00220       return rank() == cols();
00221     }
00222 
00230     inline bool isSurjective() const
00231     {
00232       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00233       return rank() == rows();
00234     }
00235 
00242     inline bool isInvertible() const
00243     {
00244       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00245       return isInjective() && isSurjective();
00246     }
00247 
00253     inline const
00254     internal::solve_retval<ColPivHouseholderQR, typename MatrixType::IdentityReturnType>
00255     inverse() const
00256     {
00257       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00258       return internal::solve_retval<ColPivHouseholderQR,typename MatrixType::IdentityReturnType>
00259                (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
00260     }
00261 
00262     inline Index rows() const { return m_qr.rows(); }
00263     inline Index cols() const { return m_qr.cols(); }
00264     const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
00265 
00283     ColPivHouseholderQR& setThreshold(const RealScalar& threshold)
00284     {
00285       m_usePrescribedThreshold = true;
00286       m_prescribedThreshold = threshold;
00287       return *this;
00288     }
00289 
00298     ColPivHouseholderQR& setThreshold(Default_t)
00299     {
00300       m_usePrescribedThreshold = false;
00301       return *this;
00302     }
00303 
00308     RealScalar threshold() const
00309     {
00310       eigen_assert(m_isInitialized || m_usePrescribedThreshold);
00311       return m_usePrescribedThreshold ? m_prescribedThreshold
00312       // this formula comes from experimenting (see "LU precision tuning" thread on the list)
00313       // and turns out to be identical to Higham's formula used already in LDLt.
00314                                       : NumTraits<Scalar>::epsilon() * m_qr.diagonalSize();
00315     }
00316 
00324     inline Index nonzeroPivots() const
00325     {
00326       eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00327       return m_nonzero_pivots;
00328     }
00329 
00333     RealScalar maxPivot() const { return m_maxpivot; }
00334 
00335   protected:
00336     MatrixType m_qr;
00337     HCoeffsType m_hCoeffs;
00338     PermutationType m_colsPermutation;
00339     IntRowVectorType m_colsTranspositions;
00340     RowVectorType m_temp;
00341     RealRowVectorType m_colSqNorms;
00342     bool m_isInitialized, m_usePrescribedThreshold;
00343     RealScalar m_prescribedThreshold, m_maxpivot;
00344     Index m_nonzero_pivots;
00345     Index m_det_pq;
00346 };
00347 
00348 template<typename MatrixType>
00349 typename MatrixType::RealScalar ColPivHouseholderQR<MatrixType>::absDeterminant() const
00350 {
00351   eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00352   eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
00353   return internal::abs(m_qr.diagonal().prod());
00354 }
00355 
00356 template<typename MatrixType>
00357 typename MatrixType::RealScalar ColPivHouseholderQR<MatrixType>::logAbsDeterminant() const
00358 {
00359   eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00360   eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
00361   return m_qr.diagonal().cwiseAbs().array().log().sum();
00362 }
00363 
00364 template<typename MatrixType>
00365 ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
00366 {
00367   Index rows = matrix.rows();
00368   Index cols = matrix.cols();
00369   Index size = matrix.diagonalSize();
00370 
00371   m_qr = matrix;
00372   m_hCoeffs.resize(size);
00373 
00374   m_temp.resize(cols);
00375 
00376   m_colsTranspositions.resize(matrix.cols());
00377   Index number_of_transpositions = 0;
00378 
00379   m_colSqNorms.resize(cols);
00380   for(Index k = 0; k < cols; ++k)
00381     m_colSqNorms.coeffRef(k) = m_qr.col(k).squaredNorm();
00382 
00383   RealScalar threshold_helper = m_colSqNorms.maxCoeff() * internal::abs2(NumTraits<Scalar>::epsilon()) / RealScalar(rows);
00384 
00385   m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
00386   m_maxpivot = RealScalar(0);
00387 
00388   for(Index k = 0; k < size; ++k)
00389   {
00390     // first, we look up in our table m_colSqNorms which column has the biggest squared norm
00391     Index biggest_col_index;
00392     RealScalar biggest_col_sq_norm = m_colSqNorms.tail(cols-k).maxCoeff(&biggest_col_index);
00393     biggest_col_index += k;
00394 
00395     // since our table m_colSqNorms accumulates imprecision at every step, we must now recompute
00396     // the actual squared norm of the selected column.
00397     // Note that not doing so does result in solve() sometimes returning inf/nan values
00398     // when running the unit test with 1000 repetitions.
00399     biggest_col_sq_norm = m_qr.col(biggest_col_index).tail(rows-k).squaredNorm();
00400 
00401     // we store that back into our table: it can't hurt to correct our table.
00402     m_colSqNorms.coeffRef(biggest_col_index) = biggest_col_sq_norm;
00403 
00404     // if the current biggest column is smaller than epsilon times the initial biggest column,
00405     // terminate to avoid generating nan/inf values.
00406     // Note that here, if we test instead for "biggest == 0", we get a failure every 1000 (or so)
00407     // repetitions of the unit test, with the result of solve() filled with large values of the order
00408     // of 1/(size*epsilon).
00409     if(biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
00410     {
00411       m_nonzero_pivots = k;
00412       m_hCoeffs.tail(size-k).setZero();
00413       m_qr.bottomRightCorner(rows-k,cols-k)
00414           .template triangularView<StrictlyLower>()
00415           .setZero();
00416       break;
00417     }
00418 
00419     // apply the transposition to the columns
00420     m_colsTranspositions.coeffRef(k) = biggest_col_index;
00421     if(k != biggest_col_index) {
00422       m_qr.col(k).swap(m_qr.col(biggest_col_index));
00423       std::swap(m_colSqNorms.coeffRef(k), m_colSqNorms.coeffRef(biggest_col_index));
00424       ++number_of_transpositions;
00425     }
00426 
00427     // generate the householder vector, store it below the diagonal
00428     RealScalar beta;
00429     m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
00430 
00431     // apply the householder transformation to the diagonal coefficient
00432     m_qr.coeffRef(k,k) = beta;
00433 
00434     // remember the maximum absolute value of diagonal coefficients
00435     if(internal::abs(beta) > m_maxpivot) m_maxpivot = internal::abs(beta);
00436 
00437     // apply the householder transformation
00438     m_qr.bottomRightCorner(rows-k, cols-k-1)
00439         .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
00440 
00441     // update our table of squared norms of the columns
00442     m_colSqNorms.tail(cols-k-1) -= m_qr.row(k).tail(cols-k-1).cwiseAbs2();
00443   }
00444 
00445   m_colsPermutation.setIdentity(PermIndexType(cols));
00446   for(PermIndexType k = 0; k < m_nonzero_pivots; ++k)
00447     m_colsPermutation.applyTranspositionOnTheRight(PermIndexType(k), PermIndexType(m_colsTranspositions.coeff(k)));
00448 
00449   m_det_pq = (number_of_transpositions%2) ? -1 : 1;
00450   m_isInitialized = true;
00451 
00452   return *this;
00453 }
00454 
00455 namespace internal {
00456 
00457 template<typename _MatrixType, typename Rhs>
00458 struct solve_retval<ColPivHouseholderQR<_MatrixType>, Rhs>
00459   : solve_retval_base<ColPivHouseholderQR<_MatrixType>, Rhs>
00460 {
00461   EIGEN_MAKE_SOLVE_HELPERS(ColPivHouseholderQR<_MatrixType>,Rhs)
00462 
00463   template<typename Dest> void evalTo(Dest& dst) const
00464   {
00465     eigen_assert(rhs().rows() == dec().rows());
00466 
00467     const int cols = dec().cols(),
00468     nonzero_pivots = dec().nonzeroPivots();
00469 
00470     if(nonzero_pivots == 0)
00471     {
00472       dst.setZero();
00473       return;
00474     }
00475 
00476     typename Rhs::PlainObject c(rhs());
00477 
00478     // Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
00479     c.applyOnTheLeft(householderSequence(dec().matrixQR(), dec().hCoeffs())
00480                      .setLength(dec().nonzeroPivots())
00481                      .transpose()
00482       );
00483 
00484     dec().matrixQR()
00485        .topLeftCorner(nonzero_pivots, nonzero_pivots)
00486        .template triangularView<Upper>()
00487        .solveInPlace(c.topRows(nonzero_pivots));
00488 
00489 
00490     typename Rhs::PlainObject d(c);
00491     d.topRows(nonzero_pivots)
00492       = dec().matrixQR()
00493        .topLeftCorner(nonzero_pivots, nonzero_pivots)
00494        .template triangularView<Upper>()
00495        * c.topRows(nonzero_pivots);
00496 
00497     for(Index i = 0; i < nonzero_pivots; ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
00498     for(Index i = nonzero_pivots; i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
00499   }
00500 };
00501 
00502 } // end namespace internal
00503 
00505 template<typename MatrixType>
00506 typename ColPivHouseholderQR<MatrixType>::HouseholderSequenceType ColPivHouseholderQR<MatrixType>
00507   ::householderQ() const
00508 {
00509   eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
00510   return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate()).setLength(m_nonzero_pivots);
00511 }
00512 
00517 template<typename Derived>
00518 const ColPivHouseholderQR<typename MatrixBase<Derived>::PlainObject>
00519 MatrixBase<Derived>::colPivHouseholderQr() const
00520 {
00521   return ColPivHouseholderQR<PlainObject>(eval());
00522 }
00523 
00524 } // end namespace Eigen
00525 
00526 #endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_H