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