Eigne Matrix Class Library
Dependents: MPC_current_control HydraulicControlBoard_SW AHRS Test_ekf ... more
src/QR/FullPivHouseholderQR.h@1:3b8049da21b8, 2019-09-24 (annotated)
- Committer:
- jsoh91
- Date:
- Tue Sep 24 00:18:23 2019 +0000
- Revision:
- 1:3b8049da21b8
- Parent:
- 0:13a5d365ba16
ignore and revise some of error parts
Who changed what in which revision?
User | Revision | Line number | New contents of line |
---|---|---|---|
ykuroda | 0:13a5d365ba16 | 1 | // This file is part of Eigen, a lightweight C++ template library |
ykuroda | 0:13a5d365ba16 | 2 | // for linear algebra. |
ykuroda | 0:13a5d365ba16 | 3 | // |
ykuroda | 0:13a5d365ba16 | 4 | // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr> |
ykuroda | 0:13a5d365ba16 | 5 | // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com> |
ykuroda | 0:13a5d365ba16 | 6 | // |
ykuroda | 0:13a5d365ba16 | 7 | // This Source Code Form is subject to the terms of the Mozilla |
ykuroda | 0:13a5d365ba16 | 8 | // Public License v. 2.0. If a copy of the MPL was not distributed |
ykuroda | 0:13a5d365ba16 | 9 | // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. |
ykuroda | 0:13a5d365ba16 | 10 | |
ykuroda | 0:13a5d365ba16 | 11 | #ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H |
ykuroda | 0:13a5d365ba16 | 12 | #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H |
ykuroda | 0:13a5d365ba16 | 13 | |
ykuroda | 0:13a5d365ba16 | 14 | namespace Eigen { |
ykuroda | 0:13a5d365ba16 | 15 | |
ykuroda | 0:13a5d365ba16 | 16 | namespace internal { |
ykuroda | 0:13a5d365ba16 | 17 | |
ykuroda | 0:13a5d365ba16 | 18 | template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType; |
ykuroda | 0:13a5d365ba16 | 19 | |
ykuroda | 0:13a5d365ba16 | 20 | template<typename MatrixType> |
ykuroda | 0:13a5d365ba16 | 21 | struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType> > |
ykuroda | 0:13a5d365ba16 | 22 | { |
ykuroda | 0:13a5d365ba16 | 23 | typedef typename MatrixType::PlainObject ReturnType; |
ykuroda | 0:13a5d365ba16 | 24 | }; |
ykuroda | 0:13a5d365ba16 | 25 | |
ykuroda | 0:13a5d365ba16 | 26 | } |
ykuroda | 0:13a5d365ba16 | 27 | |
ykuroda | 0:13a5d365ba16 | 28 | /** \ingroup QR_Module |
ykuroda | 0:13a5d365ba16 | 29 | * |
ykuroda | 0:13a5d365ba16 | 30 | * \class FullPivHouseholderQR |
ykuroda | 0:13a5d365ba16 | 31 | * |
ykuroda | 0:13a5d365ba16 | 32 | * \brief Householder rank-revealing QR decomposition of a matrix with full pivoting |
ykuroda | 0:13a5d365ba16 | 33 | * |
ykuroda | 0:13a5d365ba16 | 34 | * \param MatrixType the type of the matrix of which we are computing the QR decomposition |
ykuroda | 0:13a5d365ba16 | 35 | * |
ykuroda | 0:13a5d365ba16 | 36 | * This class performs a rank-revealing QR decomposition of a matrix \b A into matrices \b P, \b Q and \b R |
ykuroda | 0:13a5d365ba16 | 37 | * such that |
ykuroda | 0:13a5d365ba16 | 38 | * \f[ |
ykuroda | 0:13a5d365ba16 | 39 | * \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R} |
ykuroda | 0:13a5d365ba16 | 40 | * \f] |
ykuroda | 0:13a5d365ba16 | 41 | * by using Householder transformations. Here, \b P is a permutation matrix, \b Q a unitary matrix and \b R an |
ykuroda | 0:13a5d365ba16 | 42 | * upper triangular matrix. |
ykuroda | 0:13a5d365ba16 | 43 | * |
ykuroda | 0:13a5d365ba16 | 44 | * This decomposition performs a very prudent full pivoting in order to be rank-revealing and achieve optimal |
ykuroda | 0:13a5d365ba16 | 45 | * numerical stability. The trade-off is that it is slower than HouseholderQR and ColPivHouseholderQR. |
ykuroda | 0:13a5d365ba16 | 46 | * |
ykuroda | 0:13a5d365ba16 | 47 | * \sa MatrixBase::fullPivHouseholderQr() |
ykuroda | 0:13a5d365ba16 | 48 | */ |
ykuroda | 0:13a5d365ba16 | 49 | template<typename _MatrixType> class FullPivHouseholderQR |
ykuroda | 0:13a5d365ba16 | 50 | { |
ykuroda | 0:13a5d365ba16 | 51 | public: |
ykuroda | 0:13a5d365ba16 | 52 | |
ykuroda | 0:13a5d365ba16 | 53 | typedef _MatrixType MatrixType; |
ykuroda | 0:13a5d365ba16 | 54 | enum { |
ykuroda | 0:13a5d365ba16 | 55 | RowsAtCompileTime = MatrixType::RowsAtCompileTime, |
ykuroda | 0:13a5d365ba16 | 56 | ColsAtCompileTime = MatrixType::ColsAtCompileTime, |
ykuroda | 0:13a5d365ba16 | 57 | Options = MatrixType::Options, |
ykuroda | 0:13a5d365ba16 | 58 | MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, |
ykuroda | 0:13a5d365ba16 | 59 | MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime |
ykuroda | 0:13a5d365ba16 | 60 | }; |
ykuroda | 0:13a5d365ba16 | 61 | typedef typename MatrixType::Scalar Scalar; |
ykuroda | 0:13a5d365ba16 | 62 | typedef typename MatrixType::RealScalar RealScalar; |
ykuroda | 0:13a5d365ba16 | 63 | typedef typename MatrixType::Index Index; |
ykuroda | 0:13a5d365ba16 | 64 | typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType> MatrixQReturnType; |
ykuroda | 0:13a5d365ba16 | 65 | typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType; |
ykuroda | 0:13a5d365ba16 | 66 | typedef Matrix<Index, 1, |
ykuroda | 0:13a5d365ba16 | 67 | EIGEN_SIZE_MIN_PREFER_DYNAMIC(ColsAtCompileTime,RowsAtCompileTime), RowMajor, 1, |
ykuroda | 0:13a5d365ba16 | 68 | EIGEN_SIZE_MIN_PREFER_FIXED(MaxColsAtCompileTime,MaxRowsAtCompileTime)> IntDiagSizeVectorType; |
ykuroda | 0:13a5d365ba16 | 69 | typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType; |
ykuroda | 0:13a5d365ba16 | 70 | typedef typename internal::plain_row_type<MatrixType>::type RowVectorType; |
ykuroda | 0:13a5d365ba16 | 71 | typedef typename internal::plain_col_type<MatrixType>::type ColVectorType; |
ykuroda | 0:13a5d365ba16 | 72 | |
ykuroda | 0:13a5d365ba16 | 73 | /** \brief Default Constructor. |
ykuroda | 0:13a5d365ba16 | 74 | * |
ykuroda | 0:13a5d365ba16 | 75 | * The default constructor is useful in cases in which the user intends to |
ykuroda | 0:13a5d365ba16 | 76 | * perform decompositions via FullPivHouseholderQR::compute(const MatrixType&). |
ykuroda | 0:13a5d365ba16 | 77 | */ |
ykuroda | 0:13a5d365ba16 | 78 | FullPivHouseholderQR() |
ykuroda | 0:13a5d365ba16 | 79 | : m_qr(), |
ykuroda | 0:13a5d365ba16 | 80 | m_hCoeffs(), |
ykuroda | 0:13a5d365ba16 | 81 | m_rows_transpositions(), |
ykuroda | 0:13a5d365ba16 | 82 | m_cols_transpositions(), |
ykuroda | 0:13a5d365ba16 | 83 | m_cols_permutation(), |
ykuroda | 0:13a5d365ba16 | 84 | m_temp(), |
ykuroda | 0:13a5d365ba16 | 85 | m_isInitialized(false), |
ykuroda | 0:13a5d365ba16 | 86 | m_usePrescribedThreshold(false) {} |
ykuroda | 0:13a5d365ba16 | 87 | |
ykuroda | 0:13a5d365ba16 | 88 | /** \brief Default Constructor with memory preallocation |
ykuroda | 0:13a5d365ba16 | 89 | * |
ykuroda | 0:13a5d365ba16 | 90 | * Like the default constructor but with preallocation of the internal data |
ykuroda | 0:13a5d365ba16 | 91 | * according to the specified problem \a size. |
ykuroda | 0:13a5d365ba16 | 92 | * \sa FullPivHouseholderQR() |
ykuroda | 0:13a5d365ba16 | 93 | */ |
ykuroda | 0:13a5d365ba16 | 94 | FullPivHouseholderQR(Index rows, Index cols) |
ykuroda | 0:13a5d365ba16 | 95 | : m_qr(rows, cols), |
ykuroda | 0:13a5d365ba16 | 96 | m_hCoeffs((std::min)(rows,cols)), |
ykuroda | 0:13a5d365ba16 | 97 | m_rows_transpositions((std::min)(rows,cols)), |
ykuroda | 0:13a5d365ba16 | 98 | m_cols_transpositions((std::min)(rows,cols)), |
ykuroda | 0:13a5d365ba16 | 99 | m_cols_permutation(cols), |
ykuroda | 0:13a5d365ba16 | 100 | m_temp(cols), |
ykuroda | 0:13a5d365ba16 | 101 | m_isInitialized(false), |
ykuroda | 0:13a5d365ba16 | 102 | m_usePrescribedThreshold(false) {} |
ykuroda | 0:13a5d365ba16 | 103 | |
ykuroda | 0:13a5d365ba16 | 104 | /** \brief Constructs a QR factorization from a given matrix |
ykuroda | 0:13a5d365ba16 | 105 | * |
ykuroda | 0:13a5d365ba16 | 106 | * This constructor computes the QR factorization of the matrix \a matrix by calling |
ykuroda | 0:13a5d365ba16 | 107 | * the method compute(). It is a short cut for: |
ykuroda | 0:13a5d365ba16 | 108 | * |
ykuroda | 0:13a5d365ba16 | 109 | * \code |
ykuroda | 0:13a5d365ba16 | 110 | * FullPivHouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols()); |
ykuroda | 0:13a5d365ba16 | 111 | * qr.compute(matrix); |
ykuroda | 0:13a5d365ba16 | 112 | * \endcode |
ykuroda | 0:13a5d365ba16 | 113 | * |
ykuroda | 0:13a5d365ba16 | 114 | * \sa compute() |
ykuroda | 0:13a5d365ba16 | 115 | */ |
ykuroda | 0:13a5d365ba16 | 116 | FullPivHouseholderQR(const MatrixType& matrix) |
ykuroda | 0:13a5d365ba16 | 117 | : m_qr(matrix.rows(), matrix.cols()), |
ykuroda | 0:13a5d365ba16 | 118 | m_hCoeffs((std::min)(matrix.rows(), matrix.cols())), |
ykuroda | 0:13a5d365ba16 | 119 | m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())), |
ykuroda | 0:13a5d365ba16 | 120 | m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())), |
ykuroda | 0:13a5d365ba16 | 121 | m_cols_permutation(matrix.cols()), |
ykuroda | 0:13a5d365ba16 | 122 | m_temp(matrix.cols()), |
ykuroda | 0:13a5d365ba16 | 123 | m_isInitialized(false), |
ykuroda | 0:13a5d365ba16 | 124 | m_usePrescribedThreshold(false) |
ykuroda | 0:13a5d365ba16 | 125 | { |
ykuroda | 0:13a5d365ba16 | 126 | compute(matrix); |
ykuroda | 0:13a5d365ba16 | 127 | } |
ykuroda | 0:13a5d365ba16 | 128 | |
ykuroda | 0:13a5d365ba16 | 129 | /** This method finds a solution x to the equation Ax=b, where A is the matrix of which |
ykuroda | 0:13a5d365ba16 | 130 | * \c *this is the QR decomposition. |
ykuroda | 0:13a5d365ba16 | 131 | * |
ykuroda | 0:13a5d365ba16 | 132 | * \param b the right-hand-side of the equation to solve. |
ykuroda | 0:13a5d365ba16 | 133 | * |
ykuroda | 0:13a5d365ba16 | 134 | * \returns the exact or least-square solution if the rank is greater or equal to the number of columns of A, |
ykuroda | 0:13a5d365ba16 | 135 | * and an arbitrary solution otherwise. |
ykuroda | 0:13a5d365ba16 | 136 | * |
ykuroda | 0:13a5d365ba16 | 137 | * \note The case where b is a matrix is not yet implemented. Also, this |
ykuroda | 0:13a5d365ba16 | 138 | * code is space inefficient. |
ykuroda | 0:13a5d365ba16 | 139 | * |
ykuroda | 0:13a5d365ba16 | 140 | * \note_about_checking_solutions |
ykuroda | 0:13a5d365ba16 | 141 | * |
ykuroda | 0:13a5d365ba16 | 142 | * \note_about_arbitrary_choice_of_solution |
ykuroda | 0:13a5d365ba16 | 143 | * |
ykuroda | 0:13a5d365ba16 | 144 | * Example: \include FullPivHouseholderQR_solve.cpp |
ykuroda | 0:13a5d365ba16 | 145 | * Output: \verbinclude FullPivHouseholderQR_solve.out |
ykuroda | 0:13a5d365ba16 | 146 | */ |
ykuroda | 0:13a5d365ba16 | 147 | template<typename Rhs> |
ykuroda | 0:13a5d365ba16 | 148 | inline const internal::solve_retval<FullPivHouseholderQR, Rhs> |
ykuroda | 0:13a5d365ba16 | 149 | solve(const MatrixBase<Rhs>& b) const |
ykuroda | 0:13a5d365ba16 | 150 | { |
ykuroda | 0:13a5d365ba16 | 151 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 152 | return internal::solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived()); |
ykuroda | 0:13a5d365ba16 | 153 | } |
ykuroda | 0:13a5d365ba16 | 154 | |
ykuroda | 0:13a5d365ba16 | 155 | /** \returns Expression object representing the matrix Q |
ykuroda | 0:13a5d365ba16 | 156 | */ |
ykuroda | 0:13a5d365ba16 | 157 | MatrixQReturnType matrixQ(void) const; |
ykuroda | 0:13a5d365ba16 | 158 | |
ykuroda | 0:13a5d365ba16 | 159 | /** \returns a reference to the matrix where the Householder QR decomposition is stored |
ykuroda | 0:13a5d365ba16 | 160 | */ |
ykuroda | 0:13a5d365ba16 | 161 | const MatrixType& matrixQR() const |
ykuroda | 0:13a5d365ba16 | 162 | { |
ykuroda | 0:13a5d365ba16 | 163 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 164 | return m_qr; |
ykuroda | 0:13a5d365ba16 | 165 | } |
ykuroda | 0:13a5d365ba16 | 166 | |
ykuroda | 0:13a5d365ba16 | 167 | FullPivHouseholderQR& compute(const MatrixType& matrix); |
ykuroda | 0:13a5d365ba16 | 168 | |
ykuroda | 0:13a5d365ba16 | 169 | /** \returns a const reference to the column permutation matrix */ |
ykuroda | 0:13a5d365ba16 | 170 | const PermutationType& colsPermutation() const |
ykuroda | 0:13a5d365ba16 | 171 | { |
ykuroda | 0:13a5d365ba16 | 172 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 173 | return m_cols_permutation; |
ykuroda | 0:13a5d365ba16 | 174 | } |
ykuroda | 0:13a5d365ba16 | 175 | |
ykuroda | 0:13a5d365ba16 | 176 | /** \returns a const reference to the vector of indices representing the rows transpositions */ |
ykuroda | 0:13a5d365ba16 | 177 | const IntDiagSizeVectorType& rowsTranspositions() const |
ykuroda | 0:13a5d365ba16 | 178 | { |
ykuroda | 0:13a5d365ba16 | 179 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 180 | return m_rows_transpositions; |
ykuroda | 0:13a5d365ba16 | 181 | } |
ykuroda | 0:13a5d365ba16 | 182 | |
ykuroda | 0:13a5d365ba16 | 183 | /** \returns the absolute value of the determinant of the matrix of which |
ykuroda | 0:13a5d365ba16 | 184 | * *this is the QR decomposition. It has only linear complexity |
ykuroda | 0:13a5d365ba16 | 185 | * (that is, O(n) where n is the dimension of the square matrix) |
ykuroda | 0:13a5d365ba16 | 186 | * as the QR decomposition has already been computed. |
ykuroda | 0:13a5d365ba16 | 187 | * |
ykuroda | 0:13a5d365ba16 | 188 | * \note This is only for square matrices. |
ykuroda | 0:13a5d365ba16 | 189 | * |
ykuroda | 0:13a5d365ba16 | 190 | * \warning a determinant can be very big or small, so for matrices |
ykuroda | 0:13a5d365ba16 | 191 | * of large enough dimension, there is a risk of overflow/underflow. |
ykuroda | 0:13a5d365ba16 | 192 | * One way to work around that is to use logAbsDeterminant() instead. |
ykuroda | 0:13a5d365ba16 | 193 | * |
ykuroda | 0:13a5d365ba16 | 194 | * \sa logAbsDeterminant(), MatrixBase::determinant() |
ykuroda | 0:13a5d365ba16 | 195 | */ |
ykuroda | 0:13a5d365ba16 | 196 | typename MatrixType::RealScalar absDeterminant() const; |
ykuroda | 0:13a5d365ba16 | 197 | |
ykuroda | 0:13a5d365ba16 | 198 | /** \returns the natural log of the absolute value of the determinant of the matrix of which |
ykuroda | 0:13a5d365ba16 | 199 | * *this is the QR decomposition. It has only linear complexity |
ykuroda | 0:13a5d365ba16 | 200 | * (that is, O(n) where n is the dimension of the square matrix) |
ykuroda | 0:13a5d365ba16 | 201 | * as the QR decomposition has already been computed. |
ykuroda | 0:13a5d365ba16 | 202 | * |
ykuroda | 0:13a5d365ba16 | 203 | * \note This is only for square matrices. |
ykuroda | 0:13a5d365ba16 | 204 | * |
ykuroda | 0:13a5d365ba16 | 205 | * \note This method is useful to work around the risk of overflow/underflow that's inherent |
ykuroda | 0:13a5d365ba16 | 206 | * to determinant computation. |
ykuroda | 0:13a5d365ba16 | 207 | * |
ykuroda | 0:13a5d365ba16 | 208 | * \sa absDeterminant(), MatrixBase::determinant() |
ykuroda | 0:13a5d365ba16 | 209 | */ |
ykuroda | 0:13a5d365ba16 | 210 | typename MatrixType::RealScalar logAbsDeterminant() const; |
ykuroda | 0:13a5d365ba16 | 211 | |
ykuroda | 0:13a5d365ba16 | 212 | /** \returns the rank of the matrix of which *this is the QR decomposition. |
ykuroda | 0:13a5d365ba16 | 213 | * |
ykuroda | 0:13a5d365ba16 | 214 | * \note This method has to determine which pivots should be considered nonzero. |
ykuroda | 0:13a5d365ba16 | 215 | * For that, it uses the threshold value that you can control by calling |
ykuroda | 0:13a5d365ba16 | 216 | * setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 217 | */ |
ykuroda | 0:13a5d365ba16 | 218 | inline Index rank() const |
ykuroda | 0:13a5d365ba16 | 219 | { |
ykuroda | 0:13a5d365ba16 | 220 | using std::abs; |
ykuroda | 0:13a5d365ba16 | 221 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 222 | RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold(); |
ykuroda | 0:13a5d365ba16 | 223 | Index result = 0; |
ykuroda | 0:13a5d365ba16 | 224 | for(Index i = 0; i < m_nonzero_pivots; ++i) |
ykuroda | 0:13a5d365ba16 | 225 | result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold); |
ykuroda | 0:13a5d365ba16 | 226 | return result; |
ykuroda | 0:13a5d365ba16 | 227 | } |
ykuroda | 0:13a5d365ba16 | 228 | |
ykuroda | 0:13a5d365ba16 | 229 | /** \returns the dimension of the kernel of the matrix of which *this is the QR decomposition. |
ykuroda | 0:13a5d365ba16 | 230 | * |
ykuroda | 0:13a5d365ba16 | 231 | * \note This method has to determine which pivots should be considered nonzero. |
ykuroda | 0:13a5d365ba16 | 232 | * For that, it uses the threshold value that you can control by calling |
ykuroda | 0:13a5d365ba16 | 233 | * setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 234 | */ |
ykuroda | 0:13a5d365ba16 | 235 | inline Index dimensionOfKernel() const |
ykuroda | 0:13a5d365ba16 | 236 | { |
ykuroda | 0:13a5d365ba16 | 237 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 238 | return cols() - rank(); |
ykuroda | 0:13a5d365ba16 | 239 | } |
ykuroda | 0:13a5d365ba16 | 240 | |
ykuroda | 0:13a5d365ba16 | 241 | /** \returns true if the matrix of which *this is the QR decomposition represents an injective |
ykuroda | 0:13a5d365ba16 | 242 | * linear map, i.e. has trivial kernel; false otherwise. |
ykuroda | 0:13a5d365ba16 | 243 | * |
ykuroda | 0:13a5d365ba16 | 244 | * \note This method has to determine which pivots should be considered nonzero. |
ykuroda | 0:13a5d365ba16 | 245 | * For that, it uses the threshold value that you can control by calling |
ykuroda | 0:13a5d365ba16 | 246 | * setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 247 | */ |
ykuroda | 0:13a5d365ba16 | 248 | inline bool isInjective() const |
ykuroda | 0:13a5d365ba16 | 249 | { |
ykuroda | 0:13a5d365ba16 | 250 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 251 | return rank() == cols(); |
ykuroda | 0:13a5d365ba16 | 252 | } |
ykuroda | 0:13a5d365ba16 | 253 | |
ykuroda | 0:13a5d365ba16 | 254 | /** \returns true if the matrix of which *this is the QR decomposition represents a surjective |
ykuroda | 0:13a5d365ba16 | 255 | * linear map; false otherwise. |
ykuroda | 0:13a5d365ba16 | 256 | * |
ykuroda | 0:13a5d365ba16 | 257 | * \note This method has to determine which pivots should be considered nonzero. |
ykuroda | 0:13a5d365ba16 | 258 | * For that, it uses the threshold value that you can control by calling |
ykuroda | 0:13a5d365ba16 | 259 | * setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 260 | */ |
ykuroda | 0:13a5d365ba16 | 261 | inline bool isSurjective() const |
ykuroda | 0:13a5d365ba16 | 262 | { |
ykuroda | 0:13a5d365ba16 | 263 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 264 | return rank() == rows(); |
ykuroda | 0:13a5d365ba16 | 265 | } |
ykuroda | 0:13a5d365ba16 | 266 | |
ykuroda | 0:13a5d365ba16 | 267 | /** \returns true if the matrix of which *this is the QR decomposition is invertible. |
ykuroda | 0:13a5d365ba16 | 268 | * |
ykuroda | 0:13a5d365ba16 | 269 | * \note This method has to determine which pivots should be considered nonzero. |
ykuroda | 0:13a5d365ba16 | 270 | * For that, it uses the threshold value that you can control by calling |
ykuroda | 0:13a5d365ba16 | 271 | * setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 272 | */ |
ykuroda | 0:13a5d365ba16 | 273 | inline bool isInvertible() const |
ykuroda | 0:13a5d365ba16 | 274 | { |
ykuroda | 0:13a5d365ba16 | 275 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 276 | return isInjective() && isSurjective(); |
ykuroda | 0:13a5d365ba16 | 277 | } |
ykuroda | 0:13a5d365ba16 | 278 | |
ykuroda | 0:13a5d365ba16 | 279 | /** \returns the inverse of the matrix of which *this is the QR decomposition. |
ykuroda | 0:13a5d365ba16 | 280 | * |
ykuroda | 0:13a5d365ba16 | 281 | * \note If this matrix is not invertible, the returned matrix has undefined coefficients. |
ykuroda | 0:13a5d365ba16 | 282 | * Use isInvertible() to first determine whether this matrix is invertible. |
ykuroda | 0:13a5d365ba16 | 283 | */ inline const |
ykuroda | 0:13a5d365ba16 | 284 | internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType> |
ykuroda | 0:13a5d365ba16 | 285 | inverse() const |
ykuroda | 0:13a5d365ba16 | 286 | { |
ykuroda | 0:13a5d365ba16 | 287 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 288 | return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType> |
ykuroda | 0:13a5d365ba16 | 289 | (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols())); |
ykuroda | 0:13a5d365ba16 | 290 | } |
ykuroda | 0:13a5d365ba16 | 291 | |
ykuroda | 0:13a5d365ba16 | 292 | inline Index rows() const { return m_qr.rows(); } |
ykuroda | 0:13a5d365ba16 | 293 | inline Index cols() const { return m_qr.cols(); } |
ykuroda | 0:13a5d365ba16 | 294 | |
ykuroda | 0:13a5d365ba16 | 295 | /** \returns a const reference to the vector of Householder coefficients used to represent the factor \c Q. |
ykuroda | 0:13a5d365ba16 | 296 | * |
ykuroda | 0:13a5d365ba16 | 297 | * For advanced uses only. |
ykuroda | 0:13a5d365ba16 | 298 | */ |
ykuroda | 0:13a5d365ba16 | 299 | const HCoeffsType& hCoeffs() const { return m_hCoeffs; } |
ykuroda | 0:13a5d365ba16 | 300 | |
ykuroda | 0:13a5d365ba16 | 301 | /** Allows to prescribe a threshold to be used by certain methods, such as rank(), |
ykuroda | 0:13a5d365ba16 | 302 | * who need to determine when pivots are to be considered nonzero. This is not used for the |
ykuroda | 0:13a5d365ba16 | 303 | * QR decomposition itself. |
ykuroda | 0:13a5d365ba16 | 304 | * |
ykuroda | 0:13a5d365ba16 | 305 | * When it needs to get the threshold value, Eigen calls threshold(). By default, this |
ykuroda | 0:13a5d365ba16 | 306 | * uses a formula to automatically determine a reasonable threshold. |
ykuroda | 0:13a5d365ba16 | 307 | * Once you have called the present method setThreshold(const RealScalar&), |
ykuroda | 0:13a5d365ba16 | 308 | * your value is used instead. |
ykuroda | 0:13a5d365ba16 | 309 | * |
ykuroda | 0:13a5d365ba16 | 310 | * \param threshold The new value to use as the threshold. |
ykuroda | 0:13a5d365ba16 | 311 | * |
ykuroda | 0:13a5d365ba16 | 312 | * A pivot will be considered nonzero if its absolute value is strictly greater than |
ykuroda | 0:13a5d365ba16 | 313 | * \f$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \f$ |
ykuroda | 0:13a5d365ba16 | 314 | * where maxpivot is the biggest pivot. |
ykuroda | 0:13a5d365ba16 | 315 | * |
ykuroda | 0:13a5d365ba16 | 316 | * If you want to come back to the default behavior, call setThreshold(Default_t) |
ykuroda | 0:13a5d365ba16 | 317 | */ |
ykuroda | 0:13a5d365ba16 | 318 | FullPivHouseholderQR& setThreshold(const RealScalar& threshold) |
ykuroda | 0:13a5d365ba16 | 319 | { |
ykuroda | 0:13a5d365ba16 | 320 | m_usePrescribedThreshold = true; |
ykuroda | 0:13a5d365ba16 | 321 | m_prescribedThreshold = threshold; |
ykuroda | 0:13a5d365ba16 | 322 | return *this; |
ykuroda | 0:13a5d365ba16 | 323 | } |
ykuroda | 0:13a5d365ba16 | 324 | |
ykuroda | 0:13a5d365ba16 | 325 | /** Allows to come back to the default behavior, letting Eigen use its default formula for |
ykuroda | 0:13a5d365ba16 | 326 | * determining the threshold. |
ykuroda | 0:13a5d365ba16 | 327 | * |
ykuroda | 0:13a5d365ba16 | 328 | * You should pass the special object Eigen::Default as parameter here. |
ykuroda | 0:13a5d365ba16 | 329 | * \code qr.setThreshold(Eigen::Default); \endcode |
ykuroda | 0:13a5d365ba16 | 330 | * |
ykuroda | 0:13a5d365ba16 | 331 | * See the documentation of setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 332 | */ |
ykuroda | 0:13a5d365ba16 | 333 | FullPivHouseholderQR& setThreshold(Default_t) |
ykuroda | 0:13a5d365ba16 | 334 | { |
ykuroda | 0:13a5d365ba16 | 335 | m_usePrescribedThreshold = false; |
ykuroda | 0:13a5d365ba16 | 336 | return *this; |
ykuroda | 0:13a5d365ba16 | 337 | } |
ykuroda | 0:13a5d365ba16 | 338 | |
ykuroda | 0:13a5d365ba16 | 339 | /** Returns the threshold that will be used by certain methods such as rank(). |
ykuroda | 0:13a5d365ba16 | 340 | * |
ykuroda | 0:13a5d365ba16 | 341 | * See the documentation of setThreshold(const RealScalar&). |
ykuroda | 0:13a5d365ba16 | 342 | */ |
ykuroda | 0:13a5d365ba16 | 343 | RealScalar threshold() const |
ykuroda | 0:13a5d365ba16 | 344 | { |
ykuroda | 0:13a5d365ba16 | 345 | eigen_assert(m_isInitialized || m_usePrescribedThreshold); |
ykuroda | 0:13a5d365ba16 | 346 | return m_usePrescribedThreshold ? m_prescribedThreshold |
ykuroda | 0:13a5d365ba16 | 347 | // this formula comes from experimenting (see "LU precision tuning" thread on the list) |
ykuroda | 0:13a5d365ba16 | 348 | // and turns out to be identical to Higham's formula used already in LDLt. |
ykuroda | 0:13a5d365ba16 | 349 | : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize()); |
ykuroda | 0:13a5d365ba16 | 350 | } |
ykuroda | 0:13a5d365ba16 | 351 | |
ykuroda | 0:13a5d365ba16 | 352 | /** \returns the number of nonzero pivots in the QR decomposition. |
ykuroda | 0:13a5d365ba16 | 353 | * Here nonzero is meant in the exact sense, not in a fuzzy sense. |
ykuroda | 0:13a5d365ba16 | 354 | * So that notion isn't really intrinsically interesting, but it is |
ykuroda | 0:13a5d365ba16 | 355 | * still useful when implementing algorithms. |
ykuroda | 0:13a5d365ba16 | 356 | * |
ykuroda | 0:13a5d365ba16 | 357 | * \sa rank() |
ykuroda | 0:13a5d365ba16 | 358 | */ |
ykuroda | 0:13a5d365ba16 | 359 | inline Index nonzeroPivots() const |
ykuroda | 0:13a5d365ba16 | 360 | { |
ykuroda | 0:13a5d365ba16 | 361 | eigen_assert(m_isInitialized && "LU is not initialized."); |
ykuroda | 0:13a5d365ba16 | 362 | return m_nonzero_pivots; |
ykuroda | 0:13a5d365ba16 | 363 | } |
ykuroda | 0:13a5d365ba16 | 364 | |
ykuroda | 0:13a5d365ba16 | 365 | /** \returns the absolute value of the biggest pivot, i.e. the biggest |
ykuroda | 0:13a5d365ba16 | 366 | * diagonal coefficient of U. |
ykuroda | 0:13a5d365ba16 | 367 | */ |
ykuroda | 0:13a5d365ba16 | 368 | RealScalar maxPivot() const { return m_maxpivot; } |
ykuroda | 0:13a5d365ba16 | 369 | |
ykuroda | 0:13a5d365ba16 | 370 | protected: |
ykuroda | 0:13a5d365ba16 | 371 | |
ykuroda | 0:13a5d365ba16 | 372 | static void check_template_parameters() |
ykuroda | 0:13a5d365ba16 | 373 | { |
ykuroda | 0:13a5d365ba16 | 374 | EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar); |
ykuroda | 0:13a5d365ba16 | 375 | } |
ykuroda | 0:13a5d365ba16 | 376 | |
ykuroda | 0:13a5d365ba16 | 377 | MatrixType m_qr; |
ykuroda | 0:13a5d365ba16 | 378 | HCoeffsType m_hCoeffs; |
ykuroda | 0:13a5d365ba16 | 379 | IntDiagSizeVectorType m_rows_transpositions; |
ykuroda | 0:13a5d365ba16 | 380 | IntDiagSizeVectorType m_cols_transpositions; |
ykuroda | 0:13a5d365ba16 | 381 | PermutationType m_cols_permutation; |
ykuroda | 0:13a5d365ba16 | 382 | RowVectorType m_temp; |
ykuroda | 0:13a5d365ba16 | 383 | bool m_isInitialized, m_usePrescribedThreshold; |
ykuroda | 0:13a5d365ba16 | 384 | RealScalar m_prescribedThreshold, m_maxpivot; |
ykuroda | 0:13a5d365ba16 | 385 | Index m_nonzero_pivots; |
ykuroda | 0:13a5d365ba16 | 386 | RealScalar m_precision; |
ykuroda | 0:13a5d365ba16 | 387 | Index m_det_pq; |
ykuroda | 0:13a5d365ba16 | 388 | }; |
ykuroda | 0:13a5d365ba16 | 389 | |
ykuroda | 0:13a5d365ba16 | 390 | template<typename MatrixType> |
ykuroda | 0:13a5d365ba16 | 391 | typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::absDeterminant() const |
ykuroda | 0:13a5d365ba16 | 392 | { |
ykuroda | 0:13a5d365ba16 | 393 | using std::abs; |
ykuroda | 0:13a5d365ba16 | 394 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 395 | eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!"); |
ykuroda | 0:13a5d365ba16 | 396 | return abs(m_qr.diagonal().prod()); |
ykuroda | 0:13a5d365ba16 | 397 | } |
ykuroda | 0:13a5d365ba16 | 398 | |
ykuroda | 0:13a5d365ba16 | 399 | template<typename MatrixType> |
ykuroda | 0:13a5d365ba16 | 400 | typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::logAbsDeterminant() const |
ykuroda | 0:13a5d365ba16 | 401 | { |
ykuroda | 0:13a5d365ba16 | 402 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 403 | eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!"); |
ykuroda | 0:13a5d365ba16 | 404 | return m_qr.diagonal().cwiseAbs().array().log().sum(); |
ykuroda | 0:13a5d365ba16 | 405 | } |
ykuroda | 0:13a5d365ba16 | 406 | |
ykuroda | 0:13a5d365ba16 | 407 | /** Performs the QR factorization of the given matrix \a matrix. The result of |
ykuroda | 0:13a5d365ba16 | 408 | * the factorization is stored into \c *this, and a reference to \c *this |
ykuroda | 0:13a5d365ba16 | 409 | * is returned. |
ykuroda | 0:13a5d365ba16 | 410 | * |
ykuroda | 0:13a5d365ba16 | 411 | * \sa class FullPivHouseholderQR, FullPivHouseholderQR(const MatrixType&) |
ykuroda | 0:13a5d365ba16 | 412 | */ |
ykuroda | 0:13a5d365ba16 | 413 | template<typename MatrixType> |
ykuroda | 0:13a5d365ba16 | 414 | FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix) |
ykuroda | 0:13a5d365ba16 | 415 | { |
ykuroda | 0:13a5d365ba16 | 416 | check_template_parameters(); |
ykuroda | 0:13a5d365ba16 | 417 | |
ykuroda | 0:13a5d365ba16 | 418 | using std::abs; |
ykuroda | 0:13a5d365ba16 | 419 | Index rows = matrix.rows(); |
ykuroda | 0:13a5d365ba16 | 420 | Index cols = matrix.cols(); |
ykuroda | 0:13a5d365ba16 | 421 | Index size = (std::min)(rows,cols); |
ykuroda | 0:13a5d365ba16 | 422 | |
ykuroda | 0:13a5d365ba16 | 423 | m_qr = matrix; |
ykuroda | 0:13a5d365ba16 | 424 | m_hCoeffs.resize(size); |
ykuroda | 0:13a5d365ba16 | 425 | |
ykuroda | 0:13a5d365ba16 | 426 | m_temp.resize(cols); |
ykuroda | 0:13a5d365ba16 | 427 | |
ykuroda | 0:13a5d365ba16 | 428 | m_precision = NumTraits<Scalar>::epsilon() * RealScalar(size); |
ykuroda | 0:13a5d365ba16 | 429 | |
ykuroda | 0:13a5d365ba16 | 430 | m_rows_transpositions.resize(size); |
ykuroda | 0:13a5d365ba16 | 431 | m_cols_transpositions.resize(size); |
ykuroda | 0:13a5d365ba16 | 432 | Index number_of_transpositions = 0; |
ykuroda | 0:13a5d365ba16 | 433 | |
ykuroda | 0:13a5d365ba16 | 434 | RealScalar biggest(0); |
ykuroda | 0:13a5d365ba16 | 435 | |
ykuroda | 0:13a5d365ba16 | 436 | m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case) |
ykuroda | 0:13a5d365ba16 | 437 | m_maxpivot = RealScalar(0); |
ykuroda | 0:13a5d365ba16 | 438 | |
ykuroda | 0:13a5d365ba16 | 439 | for (Index k = 0; k < size; ++k) |
ykuroda | 0:13a5d365ba16 | 440 | { |
ykuroda | 0:13a5d365ba16 | 441 | Index row_of_biggest_in_corner, col_of_biggest_in_corner; |
ykuroda | 0:13a5d365ba16 | 442 | RealScalar biggest_in_corner; |
ykuroda | 0:13a5d365ba16 | 443 | |
ykuroda | 0:13a5d365ba16 | 444 | biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k) |
ykuroda | 0:13a5d365ba16 | 445 | .cwiseAbs() |
ykuroda | 0:13a5d365ba16 | 446 | .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner); |
ykuroda | 0:13a5d365ba16 | 447 | row_of_biggest_in_corner += k; |
ykuroda | 0:13a5d365ba16 | 448 | col_of_biggest_in_corner += k; |
ykuroda | 0:13a5d365ba16 | 449 | if(k==0) biggest = biggest_in_corner; |
ykuroda | 0:13a5d365ba16 | 450 | |
ykuroda | 0:13a5d365ba16 | 451 | // if the corner is negligible, then we have less than full rank, and we can finish early |
ykuroda | 0:13a5d365ba16 | 452 | if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision)) |
ykuroda | 0:13a5d365ba16 | 453 | { |
ykuroda | 0:13a5d365ba16 | 454 | m_nonzero_pivots = k; |
ykuroda | 0:13a5d365ba16 | 455 | for(Index i = k; i < size; i++) |
ykuroda | 0:13a5d365ba16 | 456 | { |
ykuroda | 0:13a5d365ba16 | 457 | m_rows_transpositions.coeffRef(i) = i; |
ykuroda | 0:13a5d365ba16 | 458 | m_cols_transpositions.coeffRef(i) = i; |
ykuroda | 0:13a5d365ba16 | 459 | m_hCoeffs.coeffRef(i) = Scalar(0); |
ykuroda | 0:13a5d365ba16 | 460 | } |
ykuroda | 0:13a5d365ba16 | 461 | break; |
ykuroda | 0:13a5d365ba16 | 462 | } |
ykuroda | 0:13a5d365ba16 | 463 | |
ykuroda | 0:13a5d365ba16 | 464 | m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner; |
ykuroda | 0:13a5d365ba16 | 465 | m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner; |
ykuroda | 0:13a5d365ba16 | 466 | if(k != row_of_biggest_in_corner) { |
ykuroda | 0:13a5d365ba16 | 467 | m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k)); |
ykuroda | 0:13a5d365ba16 | 468 | ++number_of_transpositions; |
ykuroda | 0:13a5d365ba16 | 469 | } |
ykuroda | 0:13a5d365ba16 | 470 | if(k != col_of_biggest_in_corner) { |
ykuroda | 0:13a5d365ba16 | 471 | m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner)); |
ykuroda | 0:13a5d365ba16 | 472 | ++number_of_transpositions; |
ykuroda | 0:13a5d365ba16 | 473 | } |
ykuroda | 0:13a5d365ba16 | 474 | |
ykuroda | 0:13a5d365ba16 | 475 | RealScalar beta; |
ykuroda | 0:13a5d365ba16 | 476 | m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta); |
ykuroda | 0:13a5d365ba16 | 477 | m_qr.coeffRef(k,k) = beta; |
ykuroda | 0:13a5d365ba16 | 478 | |
ykuroda | 0:13a5d365ba16 | 479 | // remember the maximum absolute value of diagonal coefficients |
ykuroda | 0:13a5d365ba16 | 480 | if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta); |
ykuroda | 0:13a5d365ba16 | 481 | |
ykuroda | 0:13a5d365ba16 | 482 | m_qr.bottomRightCorner(rows-k, cols-k-1) |
ykuroda | 0:13a5d365ba16 | 483 | .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1)); |
ykuroda | 0:13a5d365ba16 | 484 | } |
ykuroda | 0:13a5d365ba16 | 485 | |
ykuroda | 0:13a5d365ba16 | 486 | m_cols_permutation.setIdentity(cols); |
ykuroda | 0:13a5d365ba16 | 487 | for(Index k = 0; k < size; ++k) |
ykuroda | 0:13a5d365ba16 | 488 | m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k)); |
ykuroda | 0:13a5d365ba16 | 489 | |
ykuroda | 0:13a5d365ba16 | 490 | m_det_pq = (number_of_transpositions%2) ? -1 : 1; |
ykuroda | 0:13a5d365ba16 | 491 | m_isInitialized = true; |
ykuroda | 0:13a5d365ba16 | 492 | |
ykuroda | 0:13a5d365ba16 | 493 | return *this; |
ykuroda | 0:13a5d365ba16 | 494 | } |
ykuroda | 0:13a5d365ba16 | 495 | |
ykuroda | 0:13a5d365ba16 | 496 | namespace internal { |
ykuroda | 0:13a5d365ba16 | 497 | |
ykuroda | 0:13a5d365ba16 | 498 | template<typename _MatrixType, typename Rhs> |
ykuroda | 0:13a5d365ba16 | 499 | struct solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs> |
ykuroda | 0:13a5d365ba16 | 500 | : solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs> |
ykuroda | 0:13a5d365ba16 | 501 | { |
ykuroda | 0:13a5d365ba16 | 502 | EIGEN_MAKE_SOLVE_HELPERS(FullPivHouseholderQR<_MatrixType>,Rhs) |
ykuroda | 0:13a5d365ba16 | 503 | |
ykuroda | 0:13a5d365ba16 | 504 | template<typename Dest> void evalTo(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 505 | { |
ykuroda | 0:13a5d365ba16 | 506 | const Index rows = dec().rows(), cols = dec().cols(); |
ykuroda | 0:13a5d365ba16 | 507 | eigen_assert(rhs().rows() == rows); |
ykuroda | 0:13a5d365ba16 | 508 | |
ykuroda | 0:13a5d365ba16 | 509 | // FIXME introduce nonzeroPivots() and use it here. and more generally, |
ykuroda | 0:13a5d365ba16 | 510 | // make the same improvements in this dec as in FullPivLU. |
ykuroda | 0:13a5d365ba16 | 511 | if(dec().rank()==0) |
ykuroda | 0:13a5d365ba16 | 512 | { |
ykuroda | 0:13a5d365ba16 | 513 | dst.setZero(); |
ykuroda | 0:13a5d365ba16 | 514 | return; |
ykuroda | 0:13a5d365ba16 | 515 | } |
ykuroda | 0:13a5d365ba16 | 516 | |
ykuroda | 0:13a5d365ba16 | 517 | typename Rhs::PlainObject c(rhs()); |
ykuroda | 0:13a5d365ba16 | 518 | |
ykuroda | 0:13a5d365ba16 | 519 | Matrix<Scalar,1,Rhs::ColsAtCompileTime> temp(rhs().cols()); |
ykuroda | 0:13a5d365ba16 | 520 | for (Index k = 0; k < dec().rank(); ++k) |
ykuroda | 0:13a5d365ba16 | 521 | { |
ykuroda | 0:13a5d365ba16 | 522 | Index remainingSize = rows-k; |
ykuroda | 0:13a5d365ba16 | 523 | c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k))); |
ykuroda | 0:13a5d365ba16 | 524 | c.bottomRightCorner(remainingSize, rhs().cols()) |
ykuroda | 0:13a5d365ba16 | 525 | .applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1), |
ykuroda | 0:13a5d365ba16 | 526 | dec().hCoeffs().coeff(k), &temp.coeffRef(0)); |
ykuroda | 0:13a5d365ba16 | 527 | } |
ykuroda | 0:13a5d365ba16 | 528 | |
ykuroda | 0:13a5d365ba16 | 529 | dec().matrixQR() |
ykuroda | 0:13a5d365ba16 | 530 | .topLeftCorner(dec().rank(), dec().rank()) |
ykuroda | 0:13a5d365ba16 | 531 | .template triangularView<Upper>() |
ykuroda | 0:13a5d365ba16 | 532 | .solveInPlace(c.topRows(dec().rank())); |
ykuroda | 0:13a5d365ba16 | 533 | |
ykuroda | 0:13a5d365ba16 | 534 | for(Index i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i); |
ykuroda | 0:13a5d365ba16 | 535 | for(Index i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero(); |
ykuroda | 0:13a5d365ba16 | 536 | } |
ykuroda | 0:13a5d365ba16 | 537 | }; |
ykuroda | 0:13a5d365ba16 | 538 | |
ykuroda | 0:13a5d365ba16 | 539 | /** \ingroup QR_Module |
ykuroda | 0:13a5d365ba16 | 540 | * |
ykuroda | 0:13a5d365ba16 | 541 | * \brief Expression type for return value of FullPivHouseholderQR::matrixQ() |
ykuroda | 0:13a5d365ba16 | 542 | * |
ykuroda | 0:13a5d365ba16 | 543 | * \tparam MatrixType type of underlying dense matrix |
ykuroda | 0:13a5d365ba16 | 544 | */ |
ykuroda | 0:13a5d365ba16 | 545 | template<typename MatrixType> struct FullPivHouseholderQRMatrixQReturnType |
ykuroda | 0:13a5d365ba16 | 546 | : public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> > |
ykuroda | 0:13a5d365ba16 | 547 | { |
ykuroda | 0:13a5d365ba16 | 548 | public: |
ykuroda | 0:13a5d365ba16 | 549 | typedef typename MatrixType::Index Index; |
ykuroda | 0:13a5d365ba16 | 550 | typedef typename FullPivHouseholderQR<MatrixType>::IntDiagSizeVectorType IntDiagSizeVectorType; |
ykuroda | 0:13a5d365ba16 | 551 | typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType; |
ykuroda | 0:13a5d365ba16 | 552 | typedef Matrix<typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime, RowMajor, 1, |
ykuroda | 0:13a5d365ba16 | 553 | MatrixType::MaxRowsAtCompileTime> WorkVectorType; |
ykuroda | 0:13a5d365ba16 | 554 | |
ykuroda | 0:13a5d365ba16 | 555 | FullPivHouseholderQRMatrixQReturnType(const MatrixType& qr, |
ykuroda | 0:13a5d365ba16 | 556 | const HCoeffsType& hCoeffs, |
ykuroda | 0:13a5d365ba16 | 557 | const IntDiagSizeVectorType& rowsTranspositions) |
ykuroda | 0:13a5d365ba16 | 558 | : m_qr(qr), |
ykuroda | 0:13a5d365ba16 | 559 | m_hCoeffs(hCoeffs), |
ykuroda | 0:13a5d365ba16 | 560 | m_rowsTranspositions(rowsTranspositions) |
ykuroda | 0:13a5d365ba16 | 561 | {} |
ykuroda | 0:13a5d365ba16 | 562 | |
ykuroda | 0:13a5d365ba16 | 563 | template <typename ResultType> |
ykuroda | 0:13a5d365ba16 | 564 | void evalTo(ResultType& result) const |
ykuroda | 0:13a5d365ba16 | 565 | { |
ykuroda | 0:13a5d365ba16 | 566 | const Index rows = m_qr.rows(); |
ykuroda | 0:13a5d365ba16 | 567 | WorkVectorType workspace(rows); |
ykuroda | 0:13a5d365ba16 | 568 | evalTo(result, workspace); |
ykuroda | 0:13a5d365ba16 | 569 | } |
ykuroda | 0:13a5d365ba16 | 570 | |
ykuroda | 0:13a5d365ba16 | 571 | template <typename ResultType> |
ykuroda | 0:13a5d365ba16 | 572 | void evalTo(ResultType& result, WorkVectorType& workspace) const |
ykuroda | 0:13a5d365ba16 | 573 | { |
ykuroda | 0:13a5d365ba16 | 574 | using numext::conj; |
ykuroda | 0:13a5d365ba16 | 575 | // compute the product H'_0 H'_1 ... H'_n-1, |
ykuroda | 0:13a5d365ba16 | 576 | // where H_k is the k-th Householder transformation I - h_k v_k v_k' |
ykuroda | 0:13a5d365ba16 | 577 | // and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...] |
ykuroda | 0:13a5d365ba16 | 578 | const Index rows = m_qr.rows(); |
ykuroda | 0:13a5d365ba16 | 579 | const Index cols = m_qr.cols(); |
ykuroda | 0:13a5d365ba16 | 580 | const Index size = (std::min)(rows, cols); |
ykuroda | 0:13a5d365ba16 | 581 | workspace.resize(rows); |
ykuroda | 0:13a5d365ba16 | 582 | result.setIdentity(rows, rows); |
ykuroda | 0:13a5d365ba16 | 583 | for (Index k = size-1; k >= 0; k--) |
ykuroda | 0:13a5d365ba16 | 584 | { |
ykuroda | 0:13a5d365ba16 | 585 | result.block(k, k, rows-k, rows-k) |
ykuroda | 0:13a5d365ba16 | 586 | .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k)); |
ykuroda | 0:13a5d365ba16 | 587 | result.row(k).swap(result.row(m_rowsTranspositions.coeff(k))); |
ykuroda | 0:13a5d365ba16 | 588 | } |
ykuroda | 0:13a5d365ba16 | 589 | } |
ykuroda | 0:13a5d365ba16 | 590 | |
ykuroda | 0:13a5d365ba16 | 591 | Index rows() const { return m_qr.rows(); } |
ykuroda | 0:13a5d365ba16 | 592 | Index cols() const { return m_qr.rows(); } |
ykuroda | 0:13a5d365ba16 | 593 | |
ykuroda | 0:13a5d365ba16 | 594 | protected: |
ykuroda | 0:13a5d365ba16 | 595 | typename MatrixType::Nested m_qr; |
ykuroda | 0:13a5d365ba16 | 596 | typename HCoeffsType::Nested m_hCoeffs; |
ykuroda | 0:13a5d365ba16 | 597 | typename IntDiagSizeVectorType::Nested m_rowsTranspositions; |
ykuroda | 0:13a5d365ba16 | 598 | }; |
ykuroda | 0:13a5d365ba16 | 599 | |
ykuroda | 0:13a5d365ba16 | 600 | } // end namespace internal |
ykuroda | 0:13a5d365ba16 | 601 | |
ykuroda | 0:13a5d365ba16 | 602 | template<typename MatrixType> |
ykuroda | 0:13a5d365ba16 | 603 | inline typename FullPivHouseholderQR<MatrixType>::MatrixQReturnType FullPivHouseholderQR<MatrixType>::matrixQ() const |
ykuroda | 0:13a5d365ba16 | 604 | { |
ykuroda | 0:13a5d365ba16 | 605 | eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized."); |
ykuroda | 0:13a5d365ba16 | 606 | return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions); |
ykuroda | 0:13a5d365ba16 | 607 | } |
ykuroda | 0:13a5d365ba16 | 608 | |
ykuroda | 0:13a5d365ba16 | 609 | /** \return the full-pivoting Householder QR decomposition of \c *this. |
ykuroda | 0:13a5d365ba16 | 610 | * |
ykuroda | 0:13a5d365ba16 | 611 | * \sa class FullPivHouseholderQR |
ykuroda | 0:13a5d365ba16 | 612 | */ |
ykuroda | 0:13a5d365ba16 | 613 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 614 | const FullPivHouseholderQR<typename MatrixBase<Derived>::PlainObject> |
ykuroda | 0:13a5d365ba16 | 615 | MatrixBase<Derived>::fullPivHouseholderQr() const |
ykuroda | 0:13a5d365ba16 | 616 | { |
ykuroda | 0:13a5d365ba16 | 617 | return FullPivHouseholderQR<PlainObject>(eval()); |
ykuroda | 0:13a5d365ba16 | 618 | } |
ykuroda | 0:13a5d365ba16 | 619 | |
ykuroda | 0:13a5d365ba16 | 620 | } // end namespace Eigen |
ykuroda | 0:13a5d365ba16 | 621 | |
ykuroda | 0:13a5d365ba16 | 622 | #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H |