FullPivHouseholderQR.h
Go to the documentation of this file.
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_FULLPIVOTINGHOUSEHOLDERQR_H
00012 #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
00013 
00014 namespace Eigen { 
00015 
00016 namespace internal {
00017 
00018 template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType;
00019 
00020 template<typename MatrixType>
00021 struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
00022 {
00023   typedef typename MatrixType::PlainObject ReturnType;
00024 };
00025 
00026 }
00027 
00049 template<typename _MatrixType> class FullPivHouseholderQR
00050 {
00051   public:
00052 
00053     typedef _MatrixType MatrixType;
00054     enum {
00055       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00056       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00057       Options = MatrixType::Options,
00058       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00059       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00060     };
00061     typedef typename MatrixType::Scalar Scalar;
00062     typedef typename MatrixType::RealScalar RealScalar;
00063     typedef typename MatrixType::Index Index;
00064     typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType> MatrixQReturnType;
00065     typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
00066     typedef Matrix<Index, 1, ColsAtCompileTime, RowMajor, 1, MaxColsAtCompileTime> IntRowVectorType;
00067     typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
00068     typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
00069     typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
00070     typedef typename internal::plain_col_type<MatrixType>::type ColVectorType;
00071 
00077     FullPivHouseholderQR()
00078       : m_qr(),
00079         m_hCoeffs(),
00080         m_rows_transpositions(),
00081         m_cols_transpositions(),
00082         m_cols_permutation(),
00083         m_temp(),
00084         m_isInitialized(false),
00085         m_usePrescribedThreshold(false) {}
00086 
00093     FullPivHouseholderQR(Index rows, Index cols)
00094       : m_qr(rows, cols),
00095         m_hCoeffs((std::min)(rows,cols)),
00096         m_rows_transpositions(rows),
00097         m_cols_transpositions(cols),
00098         m_cols_permutation(cols),
00099         m_temp((std::min)(rows,cols)),
00100         m_isInitialized(false),
00101         m_usePrescribedThreshold(false) {}
00102 
00103     FullPivHouseholderQR(const MatrixType& matrix)
00104       : m_qr(matrix.rows(), matrix.cols()),
00105         m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
00106         m_rows_transpositions(matrix.rows()),
00107         m_cols_transpositions(matrix.cols()),
00108         m_cols_permutation(matrix.cols()),
00109         m_temp((std::min)(matrix.rows(), matrix.cols())),
00110         m_isInitialized(false),
00111         m_usePrescribedThreshold(false)
00112     {
00113       compute(matrix);
00114     }
00115 
00133     template<typename Rhs>
00134     inline const internal::solve_retval<FullPivHouseholderQR, Rhs>
00135     solve(const MatrixBase<Rhs>& b) const
00136     {
00137       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00138       return internal::solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived());
00139     }
00140 
00143     MatrixQReturnType matrixQ(void) const;
00144 
00147     const MatrixType& matrixQR() const
00148     {
00149       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00150       return m_qr;
00151     }
00152 
00153     FullPivHouseholderQR& compute(const MatrixType& matrix);
00154 
00155     const PermutationType& colsPermutation() const
00156     {
00157       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00158       return m_cols_permutation;
00159     }
00160 
00161     const IntColVectorType& rowsTranspositions() const
00162     {
00163       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00164       return m_rows_transpositions;
00165     }
00166 
00180     typename MatrixType::RealScalar absDeterminant() const;
00181 
00194     typename MatrixType::RealScalar logAbsDeterminant() const;
00195 
00202     inline Index rank() const
00203     {
00204       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00205       RealScalar premultiplied_threshold = internal::abs(m_maxpivot) * threshold();
00206       Index result = 0;
00207       for(Index i = 0; i < m_nonzero_pivots; ++i)
00208         result += (internal::abs(m_qr.coeff(i,i)) > premultiplied_threshold);
00209       return result;
00210     }
00211 
00218     inline Index dimensionOfKernel() const
00219     {
00220       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00221       return cols() - rank();
00222     }
00223 
00231     inline bool isInjective() const
00232     {
00233       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00234       return rank() == cols();
00235     }
00236 
00244     inline bool isSurjective() const
00245     {
00246       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00247       return rank() == rows();
00248     }
00249 
00256     inline bool isInvertible() const
00257     {
00258       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00259       return isInjective() && isSurjective();
00260     }
00261     inline const
00267     internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType>
00268     inverse() const
00269     {
00270       eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00271       return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType>
00272                (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
00273     }
00274 
00275     inline Index rows() const { return m_qr.rows(); }
00276     inline Index cols() const { return m_qr.cols(); }
00277     const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
00278 
00296     FullPivHouseholderQR& setThreshold(const RealScalar& threshold)
00297     {
00298       m_usePrescribedThreshold = true;
00299       m_prescribedThreshold = threshold;
00300       return *this;
00301     }
00302 
00311     FullPivHouseholderQR& setThreshold(Default_t)
00312     {
00313       m_usePrescribedThreshold = false;
00314       return *this;
00315     }
00316 
00321     RealScalar threshold() const
00322     {
00323       eigen_assert(m_isInitialized || m_usePrescribedThreshold);
00324       return m_usePrescribedThreshold ? m_prescribedThreshold
00325       // this formula comes from experimenting (see "LU precision tuning" thread on the list)
00326       // and turns out to be identical to Higham's formula used already in LDLt.
00327                                       : NumTraits<Scalar>::epsilon() * m_qr.diagonalSize();
00328     }
00329 
00337     inline Index nonzeroPivots() const
00338     {
00339       eigen_assert(m_isInitialized && "LU is not initialized.");
00340       return m_nonzero_pivots;
00341     }
00342 
00346     RealScalar maxPivot() const { return m_maxpivot; }
00347 
00348   protected:
00349     MatrixType m_qr;
00350     HCoeffsType m_hCoeffs;
00351     IntColVectorType m_rows_transpositions;
00352     IntRowVectorType m_cols_transpositions;
00353     PermutationType m_cols_permutation;
00354     RowVectorType m_temp;
00355     bool m_isInitialized, m_usePrescribedThreshold;
00356     RealScalar m_prescribedThreshold, m_maxpivot;
00357     Index m_nonzero_pivots;
00358     RealScalar m_precision;
00359     Index m_det_pq;
00360 };
00361 
00362 template<typename MatrixType>
00363 typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::absDeterminant() const
00364 {
00365   eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00366   eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
00367   return internal::abs(m_qr.diagonal().prod());
00368 }
00369 
00370 template<typename MatrixType>
00371 typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::logAbsDeterminant() const
00372 {
00373   eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00374   eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
00375   return m_qr.diagonal().cwiseAbs().array().log().sum();
00376 }
00377 
00378 template<typename MatrixType>
00379 FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
00380 {
00381   Index rows = matrix.rows();
00382   Index cols = matrix.cols();
00383   Index size = (std::min)(rows,cols);
00384 
00385   m_qr = matrix;
00386   m_hCoeffs.resize(size);
00387 
00388   m_temp.resize(cols);
00389 
00390   m_precision = NumTraits<Scalar>::epsilon() * size;
00391 
00392   m_rows_transpositions.resize(matrix.rows());
00393   m_cols_transpositions.resize(matrix.cols());
00394   Index number_of_transpositions = 0;
00395 
00396   RealScalar biggest(0);
00397 
00398   m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
00399   m_maxpivot = RealScalar(0);
00400 
00401   for (Index k = 0; k < size; ++k)
00402   {
00403     Index row_of_biggest_in_corner, col_of_biggest_in_corner;
00404     RealScalar biggest_in_corner;
00405 
00406     biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
00407                             .cwiseAbs()
00408                             .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
00409     row_of_biggest_in_corner += k;
00410     col_of_biggest_in_corner += k;
00411     if(k==0) biggest = biggest_in_corner;
00412 
00413     // if the corner is negligible, then we have less than full rank, and we can finish early
00414     if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
00415     {
00416       m_nonzero_pivots = k;
00417       for(Index i = k; i < size; i++)
00418       {
00419         m_rows_transpositions.coeffRef(i) = i;
00420         m_cols_transpositions.coeffRef(i) = i;
00421         m_hCoeffs.coeffRef(i) = Scalar(0);
00422       }
00423       break;
00424     }
00425 
00426     m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
00427     m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
00428     if(k != row_of_biggest_in_corner) {
00429       m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
00430       ++number_of_transpositions;
00431     }
00432     if(k != col_of_biggest_in_corner) {
00433       m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
00434       ++number_of_transpositions;
00435     }
00436 
00437     RealScalar beta;
00438     m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
00439     m_qr.coeffRef(k,k) = beta;
00440 
00441     // remember the maximum absolute value of diagonal coefficients
00442     if(internal::abs(beta) > m_maxpivot) m_maxpivot = internal::abs(beta);
00443 
00444     m_qr.bottomRightCorner(rows-k, cols-k-1)
00445         .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
00446   }
00447 
00448   m_cols_permutation.setIdentity(cols);
00449   for(Index k = 0; k < size; ++k)
00450     m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
00451 
00452   m_det_pq = (number_of_transpositions%2) ? -1 : 1;
00453   m_isInitialized = true;
00454 
00455   return *this;
00456 }
00457 
00458 namespace internal {
00459 
00460 template<typename _MatrixType, typename Rhs>
00461 struct solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
00462   : solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs>
00463 {
00464   EIGEN_MAKE_SOLVE_HELPERS(FullPivHouseholderQR<_MatrixType>,Rhs)
00465 
00466   template<typename Dest> void evalTo(Dest& dst) const
00467   {
00468     const Index rows = dec().rows(), cols = dec().cols();
00469     eigen_assert(rhs().rows() == rows);
00470 
00471     // FIXME introduce nonzeroPivots() and use it here. and more generally,
00472     // make the same improvements in this dec as in FullPivLU.
00473     if(dec().rank()==0)
00474     {
00475       dst.setZero();
00476       return;
00477     }
00478 
00479     typename Rhs::PlainObject c(rhs());
00480 
00481     Matrix<Scalar,1,Rhs::ColsAtCompileTime> temp(rhs().cols());
00482     for (Index k = 0; k < dec().rank(); ++k)
00483     {
00484       Index remainingSize = rows-k;
00485       c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k)));
00486       c.bottomRightCorner(remainingSize, rhs().cols())
00487        .applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1),
00488                                   dec().hCoeffs().coeff(k), &temp.coeffRef(0));
00489     }
00490 
00491     if(!dec().isSurjective())
00492     {
00493       // is c is in the image of R ?
00494       RealScalar biggest_in_upper_part_of_c = c.topRows(   dec().rank()     ).cwiseAbs().maxCoeff();
00495       RealScalar biggest_in_lower_part_of_c = c.bottomRows(rows-dec().rank()).cwiseAbs().maxCoeff();
00496       // FIXME brain dead
00497       const RealScalar m_precision = NumTraits<Scalar>::epsilon() * (std::min)(rows,cols);
00498       // this internal:: prefix is needed by at least gcc 3.4 and ICC
00499       if(!internal::isMuchSmallerThan(biggest_in_lower_part_of_c, biggest_in_upper_part_of_c, m_precision))
00500         return;
00501     }
00502     dec().matrixQR()
00503        .topLeftCorner(dec().rank(), dec().rank())
00504        .template triangularView<Upper>()
00505        .solveInPlace(c.topRows(dec().rank()));
00506 
00507     for(Index i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
00508     for(Index i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
00509   }
00510 };
00511 
00518 template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType
00519   : public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
00520 {
00521 public:
00522   typedef typename MatrixType::Index Index;
00523   typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
00524   typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
00525   typedef Matrix<typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime, RowMajor, 1,
00526                  MatrixType::MaxRowsAtCompileTime> WorkVectorType;
00527 
00528   FullPivHouseholderQRMatrixQReturnType(const MatrixType&       qr,
00529                                         const HCoeffsType&      hCoeffs,
00530                                         const IntColVectorType& rowsTranspositions)
00531     : m_qr(qr),
00532       m_hCoeffs(hCoeffs),
00533       m_rowsTranspositions(rowsTranspositions)
00534       {}
00535 
00536   template <typename ResultType>
00537   void evalTo(ResultType& result) const
00538   {
00539     const Index rows = m_qr.rows();
00540     WorkVectorType workspace(rows);
00541     evalTo(result, workspace);
00542   }
00543 
00544   template <typename ResultType>
00545   void evalTo(ResultType& result, WorkVectorType& workspace) const
00546   {
00547     // compute the product H'_0 H'_1 ... H'_n-1,
00548     // where H_k is the k-th Householder transformation I - h_k v_k v_k'
00549     // and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
00550     const Index rows = m_qr.rows();
00551     const Index cols = m_qr.cols();
00552     const Index size = (std::min)(rows, cols);
00553     workspace.resize(rows);
00554     result.setIdentity(rows, rows);
00555     for (Index k = size-1; k >= 0; k--)
00556     {
00557       result.block(k, k, rows-k, rows-k)
00558             .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), internal::conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k));
00559       result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
00560     }
00561   }
00562 
00563     Index rows() const { return m_qr.rows(); }
00564     Index cols() const { return m_qr.rows(); }
00565 
00566 protected:
00567   typename MatrixType::Nested m_qr;
00568   typename HCoeffsType::Nested m_hCoeffs;
00569   typename IntColVectorType::Nested m_rowsTranspositions;
00570 };
00571 
00572 } // end namespace internal
00573 
00574 template<typename MatrixType>
00575 inline typename FullPivHouseholderQR<MatrixType>::MatrixQReturnType FullPivHouseholderQR<MatrixType>::matrixQ() const
00576 {
00577   eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
00578   return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
00579 }
00580 
00585 template<typename Derived>
00586 const FullPivHouseholderQR<typename MatrixBase<Derived>::PlainObject>
00587 MatrixBase<Derived>::fullPivHouseholderQr() const
00588 {
00589   return FullPivHouseholderQR<PlainObject>(eval());
00590 }
00591 
00592 } // end namespace Eigen
00593 
00594 #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Mon Oct 6 2014 12:24:29