CompleteOrthogonalDecomposition.h 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2016 Rasmus Munk Larsen <rmlarsen@google.com>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
  10. #define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
  11. namespace Eigen {
  12. namespace internal {
  13. template <typename _MatrixType>
  14. struct traits<CompleteOrthogonalDecomposition<_MatrixType> >
  15. : traits<_MatrixType> {
  16. typedef MatrixXpr XprKind;
  17. typedef SolverStorage StorageKind;
  18. typedef int StorageIndex;
  19. enum { Flags = 0 };
  20. };
  21. } // end namespace internal
  22. /** \ingroup QR_Module
  23. *
  24. * \class CompleteOrthogonalDecomposition
  25. *
  26. * \brief Complete orthogonal decomposition (COD) of a matrix.
  27. *
  28. * \param MatrixType the type of the matrix of which we are computing the COD.
  29. *
  30. * This class performs a rank-revealing complete orthogonal decomposition of a
  31. * matrix \b A into matrices \b P, \b Q, \b T, and \b Z such that
  32. * \f[
  33. * \mathbf{A} \, \mathbf{P} = \mathbf{Q} \,
  34. * \begin{bmatrix} \mathbf{T} & \mathbf{0} \\
  35. * \mathbf{0} & \mathbf{0} \end{bmatrix} \, \mathbf{Z}
  36. * \f]
  37. * by using Householder transformations. Here, \b P is a permutation matrix,
  38. * \b Q and \b Z are unitary matrices and \b T an upper triangular matrix of
  39. * size rank-by-rank. \b A may be rank deficient.
  40. *
  41. * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
  42. *
  43. * \sa MatrixBase::completeOrthogonalDecomposition()
  44. */
  45. template <typename _MatrixType> class CompleteOrthogonalDecomposition
  46. : public SolverBase<CompleteOrthogonalDecomposition<_MatrixType> >
  47. {
  48. public:
  49. typedef _MatrixType MatrixType;
  50. typedef SolverBase<CompleteOrthogonalDecomposition> Base;
  51. template<typename Derived>
  52. friend struct internal::solve_assertion;
  53. EIGEN_GENERIC_PUBLIC_INTERFACE(CompleteOrthogonalDecomposition)
  54. enum {
  55. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  56. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  57. };
  58. typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
  59. typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime>
  60. PermutationType;
  61. typedef typename internal::plain_row_type<MatrixType, Index>::type
  62. IntRowVectorType;
  63. typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
  64. typedef typename internal::plain_row_type<MatrixType, RealScalar>::type
  65. RealRowVectorType;
  66. typedef HouseholderSequence<
  67. MatrixType, typename internal::remove_all<
  68. typename HCoeffsType::ConjugateReturnType>::type>
  69. HouseholderSequenceType;
  70. typedef typename MatrixType::PlainObject PlainObject;
  71. private:
  72. typedef typename PermutationType::Index PermIndexType;
  73. public:
  74. /**
  75. * \brief Default Constructor.
  76. *
  77. * The default constructor is useful in cases in which the user intends to
  78. * perform decompositions via
  79. * \c CompleteOrthogonalDecomposition::compute(const* MatrixType&).
  80. */
  81. CompleteOrthogonalDecomposition() : m_cpqr(), m_zCoeffs(), m_temp() {}
  82. /** \brief Default Constructor with memory preallocation
  83. *
  84. * Like the default constructor but with preallocation of the internal data
  85. * according to the specified problem \a size.
  86. * \sa CompleteOrthogonalDecomposition()
  87. */
  88. CompleteOrthogonalDecomposition(Index rows, Index cols)
  89. : m_cpqr(rows, cols), m_zCoeffs((std::min)(rows, cols)), m_temp(cols) {}
  90. /** \brief Constructs a complete orthogonal decomposition from a given
  91. * matrix.
  92. *
  93. * This constructor computes the complete orthogonal decomposition of the
  94. * matrix \a matrix by calling the method compute(). The default
  95. * threshold for rank determination will be used. It is a short cut for:
  96. *
  97. * \code
  98. * CompleteOrthogonalDecomposition<MatrixType> cod(matrix.rows(),
  99. * matrix.cols());
  100. * cod.setThreshold(Default);
  101. * cod.compute(matrix);
  102. * \endcode
  103. *
  104. * \sa compute()
  105. */
  106. template <typename InputType>
  107. explicit CompleteOrthogonalDecomposition(const EigenBase<InputType>& matrix)
  108. : m_cpqr(matrix.rows(), matrix.cols()),
  109. m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
  110. m_temp(matrix.cols())
  111. {
  112. compute(matrix.derived());
  113. }
  114. /** \brief Constructs a complete orthogonal decomposition from a given matrix
  115. *
  116. * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref.
  117. *
  118. * \sa CompleteOrthogonalDecomposition(const EigenBase&)
  119. */
  120. template<typename InputType>
  121. explicit CompleteOrthogonalDecomposition(EigenBase<InputType>& matrix)
  122. : m_cpqr(matrix.derived()),
  123. m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
  124. m_temp(matrix.cols())
  125. {
  126. computeInPlace();
  127. }
  128. #ifdef EIGEN_PARSED_BY_DOXYGEN
  129. /** This method computes the minimum-norm solution X to a least squares
  130. * problem \f[\mathrm{minimize} \|A X - B\|, \f] where \b A is the matrix of
  131. * which \c *this is the complete orthogonal decomposition.
  132. *
  133. * \param b the right-hand sides of the problem to solve.
  134. *
  135. * \returns a solution.
  136. *
  137. */
  138. template <typename Rhs>
  139. inline const Solve<CompleteOrthogonalDecomposition, Rhs> solve(
  140. const MatrixBase<Rhs>& b) const;
  141. #endif
  142. HouseholderSequenceType householderQ(void) const;
  143. HouseholderSequenceType matrixQ(void) const { return m_cpqr.householderQ(); }
  144. /** \returns the matrix \b Z.
  145. */
  146. MatrixType matrixZ() const {
  147. MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
  148. applyZOnTheLeftInPlace<false>(Z);
  149. return Z;
  150. }
  151. /** \returns a reference to the matrix where the complete orthogonal
  152. * decomposition is stored
  153. */
  154. const MatrixType& matrixQTZ() const { return m_cpqr.matrixQR(); }
  155. /** \returns a reference to the matrix where the complete orthogonal
  156. * decomposition is stored.
  157. * \warning The strict lower part and \code cols() - rank() \endcode right
  158. * columns of this matrix contains internal values.
  159. * Only the upper triangular part should be referenced. To get it, use
  160. * \code matrixT().template triangularView<Upper>() \endcode
  161. * For rank-deficient matrices, use
  162. * \code
  163. * matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()
  164. * \endcode
  165. */
  166. const MatrixType& matrixT() const { return m_cpqr.matrixQR(); }
  167. template <typename InputType>
  168. CompleteOrthogonalDecomposition& compute(const EigenBase<InputType>& matrix) {
  169. // Compute the column pivoted QR factorization A P = Q R.
  170. m_cpqr.compute(matrix);
  171. computeInPlace();
  172. return *this;
  173. }
  174. /** \returns a const reference to the column permutation matrix */
  175. const PermutationType& colsPermutation() const {
  176. return m_cpqr.colsPermutation();
  177. }
  178. /** \returns the absolute value of the determinant of the matrix of which
  179. * *this is the complete orthogonal decomposition. It has only linear
  180. * complexity (that is, O(n) where n is the dimension of the square matrix)
  181. * as the complete orthogonal decomposition has already been computed.
  182. *
  183. * \note This is only for square matrices.
  184. *
  185. * \warning a determinant can be very big or small, so for matrices
  186. * of large enough dimension, there is a risk of overflow/underflow.
  187. * One way to work around that is to use logAbsDeterminant() instead.
  188. *
  189. * \sa logAbsDeterminant(), MatrixBase::determinant()
  190. */
  191. typename MatrixType::RealScalar absDeterminant() const;
  192. /** \returns the natural log of the absolute value of the determinant of the
  193. * matrix of which *this is the complete orthogonal decomposition. It has
  194. * only linear complexity (that is, O(n) where n is the dimension of the
  195. * square matrix) as the complete orthogonal decomposition has already been
  196. * computed.
  197. *
  198. * \note This is only for square matrices.
  199. *
  200. * \note This method is useful to work around the risk of overflow/underflow
  201. * that's inherent to determinant computation.
  202. *
  203. * \sa absDeterminant(), MatrixBase::determinant()
  204. */
  205. typename MatrixType::RealScalar logAbsDeterminant() const;
  206. /** \returns the rank of the matrix of which *this is the complete orthogonal
  207. * decomposition.
  208. *
  209. * \note This method has to determine which pivots should be considered
  210. * nonzero. For that, it uses the threshold value that you can control by
  211. * calling setThreshold(const RealScalar&).
  212. */
  213. inline Index rank() const { return m_cpqr.rank(); }
  214. /** \returns the dimension of the kernel of the matrix of which *this is the
  215. * complete orthogonal decomposition.
  216. *
  217. * \note This method has to determine which pivots should be considered
  218. * nonzero. For that, it uses the threshold value that you can control by
  219. * calling setThreshold(const RealScalar&).
  220. */
  221. inline Index dimensionOfKernel() const { return m_cpqr.dimensionOfKernel(); }
  222. /** \returns true if the matrix of which *this is the decomposition represents
  223. * an injective linear map, i.e. has trivial kernel; false otherwise.
  224. *
  225. * \note This method has to determine which pivots should be considered
  226. * nonzero. For that, it uses the threshold value that you can control by
  227. * calling setThreshold(const RealScalar&).
  228. */
  229. inline bool isInjective() const { return m_cpqr.isInjective(); }
  230. /** \returns true if the matrix of which *this is the decomposition represents
  231. * a surjective linear map; false otherwise.
  232. *
  233. * \note This method has to determine which pivots should be considered
  234. * nonzero. For that, it uses the threshold value that you can control by
  235. * calling setThreshold(const RealScalar&).
  236. */
  237. inline bool isSurjective() const { return m_cpqr.isSurjective(); }
  238. /** \returns true if the matrix of which *this is the complete orthogonal
  239. * decomposition is invertible.
  240. *
  241. * \note This method has to determine which pivots should be considered
  242. * nonzero. For that, it uses the threshold value that you can control by
  243. * calling setThreshold(const RealScalar&).
  244. */
  245. inline bool isInvertible() const { return m_cpqr.isInvertible(); }
  246. /** \returns the pseudo-inverse of the matrix of which *this is the complete
  247. * orthogonal decomposition.
  248. * \warning: Do not compute \c this->pseudoInverse()*rhs to solve a linear systems.
  249. * It is more efficient and numerically stable to call \c this->solve(rhs).
  250. */
  251. inline const Inverse<CompleteOrthogonalDecomposition> pseudoInverse() const
  252. {
  253. eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
  254. return Inverse<CompleteOrthogonalDecomposition>(*this);
  255. }
  256. inline Index rows() const { return m_cpqr.rows(); }
  257. inline Index cols() const { return m_cpqr.cols(); }
  258. /** \returns a const reference to the vector of Householder coefficients used
  259. * to represent the factor \c Q.
  260. *
  261. * For advanced uses only.
  262. */
  263. inline const HCoeffsType& hCoeffs() const { return m_cpqr.hCoeffs(); }
  264. /** \returns a const reference to the vector of Householder coefficients
  265. * used to represent the factor \c Z.
  266. *
  267. * For advanced uses only.
  268. */
  269. const HCoeffsType& zCoeffs() const { return m_zCoeffs; }
  270. /** Allows to prescribe a threshold to be used by certain methods, such as
  271. * rank(), who need to determine when pivots are to be considered nonzero.
  272. * Most be called before calling compute().
  273. *
  274. * When it needs to get the threshold value, Eigen calls threshold(). By
  275. * default, this uses a formula to automatically determine a reasonable
  276. * threshold. Once you have called the present method
  277. * setThreshold(const RealScalar&), your value is used instead.
  278. *
  279. * \param threshold The new value to use as the threshold.
  280. *
  281. * A pivot will be considered nonzero if its absolute value is strictly
  282. * greater than
  283. * \f$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \f$
  284. * where maxpivot is the biggest pivot.
  285. *
  286. * If you want to come back to the default behavior, call
  287. * setThreshold(Default_t)
  288. */
  289. CompleteOrthogonalDecomposition& setThreshold(const RealScalar& threshold) {
  290. m_cpqr.setThreshold(threshold);
  291. return *this;
  292. }
  293. /** Allows to come back to the default behavior, letting Eigen use its default
  294. * formula for determining the threshold.
  295. *
  296. * You should pass the special object Eigen::Default as parameter here.
  297. * \code qr.setThreshold(Eigen::Default); \endcode
  298. *
  299. * See the documentation of setThreshold(const RealScalar&).
  300. */
  301. CompleteOrthogonalDecomposition& setThreshold(Default_t) {
  302. m_cpqr.setThreshold(Default);
  303. return *this;
  304. }
  305. /** Returns the threshold that will be used by certain methods such as rank().
  306. *
  307. * See the documentation of setThreshold(const RealScalar&).
  308. */
  309. RealScalar threshold() const { return m_cpqr.threshold(); }
  310. /** \returns the number of nonzero pivots in the complete orthogonal
  311. * decomposition. Here nonzero is meant in the exact sense, not in a
  312. * fuzzy sense. So that notion isn't really intrinsically interesting,
  313. * but it is still useful when implementing algorithms.
  314. *
  315. * \sa rank()
  316. */
  317. inline Index nonzeroPivots() const { return m_cpqr.nonzeroPivots(); }
  318. /** \returns the absolute value of the biggest pivot, i.e. the biggest
  319. * diagonal coefficient of R.
  320. */
  321. inline RealScalar maxPivot() const { return m_cpqr.maxPivot(); }
  322. /** \brief Reports whether the complete orthogonal decomposition was
  323. * successful.
  324. *
  325. * \note This function always returns \c Success. It is provided for
  326. * compatibility
  327. * with other factorization routines.
  328. * \returns \c Success
  329. */
  330. ComputationInfo info() const {
  331. eigen_assert(m_cpqr.m_isInitialized && "Decomposition is not initialized.");
  332. return Success;
  333. }
  334. #ifndef EIGEN_PARSED_BY_DOXYGEN
  335. template <typename RhsType, typename DstType>
  336. void _solve_impl(const RhsType& rhs, DstType& dst) const;
  337. template<bool Conjugate, typename RhsType, typename DstType>
  338. void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
  339. #endif
  340. protected:
  341. static void check_template_parameters() {
  342. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  343. }
  344. template<bool Transpose_, typename Rhs>
  345. void _check_solve_assertion(const Rhs& b) const {
  346. EIGEN_ONLY_USED_FOR_DEBUG(b);
  347. eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
  348. eigen_assert((Transpose_?derived().cols():derived().rows())==b.rows() && "CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
  349. }
  350. void computeInPlace();
  351. /** Overwrites \b rhs with \f$ \mathbf{Z} * \mathbf{rhs} \f$ or
  352. * \f$ \mathbf{\overline Z} * \mathbf{rhs} \f$ if \c Conjugate
  353. * is set to \c true.
  354. */
  355. template <bool Conjugate, typename Rhs>
  356. void applyZOnTheLeftInPlace(Rhs& rhs) const;
  357. /** Overwrites \b rhs with \f$ \mathbf{Z}^* * \mathbf{rhs} \f$.
  358. */
  359. template <typename Rhs>
  360. void applyZAdjointOnTheLeftInPlace(Rhs& rhs) const;
  361. ColPivHouseholderQR<MatrixType> m_cpqr;
  362. HCoeffsType m_zCoeffs;
  363. RowVectorType m_temp;
  364. };
  365. template <typename MatrixType>
  366. typename MatrixType::RealScalar
  367. CompleteOrthogonalDecomposition<MatrixType>::absDeterminant() const {
  368. return m_cpqr.absDeterminant();
  369. }
  370. template <typename MatrixType>
  371. typename MatrixType::RealScalar
  372. CompleteOrthogonalDecomposition<MatrixType>::logAbsDeterminant() const {
  373. return m_cpqr.logAbsDeterminant();
  374. }
  375. /** Performs the complete orthogonal decomposition of the given matrix \a
  376. * matrix. The result of the factorization is stored into \c *this, and a
  377. * reference to \c *this is returned.
  378. *
  379. * \sa class CompleteOrthogonalDecomposition,
  380. * CompleteOrthogonalDecomposition(const MatrixType&)
  381. */
  382. template <typename MatrixType>
  383. void CompleteOrthogonalDecomposition<MatrixType>::computeInPlace()
  384. {
  385. check_template_parameters();
  386. // the column permutation is stored as int indices, so just to be sure:
  387. eigen_assert(m_cpqr.cols() <= NumTraits<int>::highest());
  388. const Index rank = m_cpqr.rank();
  389. const Index cols = m_cpqr.cols();
  390. const Index rows = m_cpqr.rows();
  391. m_zCoeffs.resize((std::min)(rows, cols));
  392. m_temp.resize(cols);
  393. if (rank < cols) {
  394. // We have reduced the (permuted) matrix to the form
  395. // [R11 R12]
  396. // [ 0 R22]
  397. // where R11 is r-by-r (r = rank) upper triangular, R12 is
  398. // r-by-(n-r), and R22 is empty or the norm of R22 is negligible.
  399. // We now compute the complete orthogonal decomposition by applying
  400. // Householder transformations from the right to the upper trapezoidal
  401. // matrix X = [R11 R12] to zero out R12 and obtain the factorization
  402. // [R11 R12] = [T11 0] * Z, where T11 is r-by-r upper triangular and
  403. // Z = Z(0) * Z(1) ... Z(r-1) is an n-by-n orthogonal matrix.
  404. // We store the data representing Z in R12 and m_zCoeffs.
  405. for (Index k = rank - 1; k >= 0; --k) {
  406. if (k != rank - 1) {
  407. // Given the API for Householder reflectors, it is more convenient if
  408. // we swap the leading parts of columns k and r-1 (zero-based) to form
  409. // the matrix X_k = [X(0:k, k), X(0:k, r:n)]
  410. m_cpqr.m_qr.col(k).head(k + 1).swap(
  411. m_cpqr.m_qr.col(rank - 1).head(k + 1));
  412. }
  413. // Construct Householder reflector Z(k) to zero out the last row of X_k,
  414. // i.e. choose Z(k) such that
  415. // [X(k, k), X(k, r:n)] * Z(k) = [beta, 0, .., 0].
  416. RealScalar beta;
  417. m_cpqr.m_qr.row(k)
  418. .tail(cols - rank + 1)
  419. .makeHouseholderInPlace(m_zCoeffs(k), beta);
  420. m_cpqr.m_qr(k, rank - 1) = beta;
  421. if (k > 0) {
  422. // Apply Z(k) to the first k rows of X_k
  423. m_cpqr.m_qr.topRightCorner(k, cols - rank + 1)
  424. .applyHouseholderOnTheRight(
  425. m_cpqr.m_qr.row(k).tail(cols - rank).adjoint(), m_zCoeffs(k),
  426. &m_temp(0));
  427. }
  428. if (k != rank - 1) {
  429. // Swap X(0:k,k) back to its proper location.
  430. m_cpqr.m_qr.col(k).head(k + 1).swap(
  431. m_cpqr.m_qr.col(rank - 1).head(k + 1));
  432. }
  433. }
  434. }
  435. }
  436. template <typename MatrixType>
  437. template <bool Conjugate, typename Rhs>
  438. void CompleteOrthogonalDecomposition<MatrixType>::applyZOnTheLeftInPlace(
  439. Rhs& rhs) const {
  440. const Index cols = this->cols();
  441. const Index nrhs = rhs.cols();
  442. const Index rank = this->rank();
  443. Matrix<typename Rhs::Scalar, Dynamic, 1> temp((std::max)(cols, nrhs));
  444. for (Index k = rank-1; k >= 0; --k) {
  445. if (k != rank - 1) {
  446. rhs.row(k).swap(rhs.row(rank - 1));
  447. }
  448. rhs.middleRows(rank - 1, cols - rank + 1)
  449. .applyHouseholderOnTheLeft(
  450. matrixQTZ().row(k).tail(cols - rank).transpose().template conjugateIf<!Conjugate>(), zCoeffs().template conjugateIf<Conjugate>()(k),
  451. &temp(0));
  452. if (k != rank - 1) {
  453. rhs.row(k).swap(rhs.row(rank - 1));
  454. }
  455. }
  456. }
  457. template <typename MatrixType>
  458. template <typename Rhs>
  459. void CompleteOrthogonalDecomposition<MatrixType>::applyZAdjointOnTheLeftInPlace(
  460. Rhs& rhs) const {
  461. const Index cols = this->cols();
  462. const Index nrhs = rhs.cols();
  463. const Index rank = this->rank();
  464. Matrix<typename Rhs::Scalar, Dynamic, 1> temp((std::max)(cols, nrhs));
  465. for (Index k = 0; k < rank; ++k) {
  466. if (k != rank - 1) {
  467. rhs.row(k).swap(rhs.row(rank - 1));
  468. }
  469. rhs.middleRows(rank - 1, cols - rank + 1)
  470. .applyHouseholderOnTheLeft(
  471. matrixQTZ().row(k).tail(cols - rank).adjoint(), zCoeffs()(k),
  472. &temp(0));
  473. if (k != rank - 1) {
  474. rhs.row(k).swap(rhs.row(rank - 1));
  475. }
  476. }
  477. }
  478. #ifndef EIGEN_PARSED_BY_DOXYGEN
  479. template <typename _MatrixType>
  480. template <typename RhsType, typename DstType>
  481. void CompleteOrthogonalDecomposition<_MatrixType>::_solve_impl(
  482. const RhsType& rhs, DstType& dst) const {
  483. const Index rank = this->rank();
  484. if (rank == 0) {
  485. dst.setZero();
  486. return;
  487. }
  488. // Compute c = Q^* * rhs
  489. typename RhsType::PlainObject c(rhs);
  490. c.applyOnTheLeft(matrixQ().setLength(rank).adjoint());
  491. // Solve T z = c(1:rank, :)
  492. dst.topRows(rank) = matrixT()
  493. .topLeftCorner(rank, rank)
  494. .template triangularView<Upper>()
  495. .solve(c.topRows(rank));
  496. const Index cols = this->cols();
  497. if (rank < cols) {
  498. // Compute y = Z^* * [ z ]
  499. // [ 0 ]
  500. dst.bottomRows(cols - rank).setZero();
  501. applyZAdjointOnTheLeftInPlace(dst);
  502. }
  503. // Undo permutation to get x = P^{-1} * y.
  504. dst = colsPermutation() * dst;
  505. }
  506. template<typename _MatrixType>
  507. template<bool Conjugate, typename RhsType, typename DstType>
  508. void CompleteOrthogonalDecomposition<_MatrixType>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
  509. {
  510. const Index rank = this->rank();
  511. if (rank == 0) {
  512. dst.setZero();
  513. return;
  514. }
  515. typename RhsType::PlainObject c(colsPermutation().transpose()*rhs);
  516. if (rank < cols()) {
  517. applyZOnTheLeftInPlace<!Conjugate>(c);
  518. }
  519. matrixT().topLeftCorner(rank, rank)
  520. .template triangularView<Upper>()
  521. .transpose().template conjugateIf<Conjugate>()
  522. .solveInPlace(c.topRows(rank));
  523. dst.topRows(rank) = c.topRows(rank);
  524. dst.bottomRows(rows()-rank).setZero();
  525. dst.applyOnTheLeft(householderQ().setLength(rank).template conjugateIf<!Conjugate>() );
  526. }
  527. #endif
  528. namespace internal {
  529. template<typename MatrixType>
  530. struct traits<Inverse<CompleteOrthogonalDecomposition<MatrixType> > >
  531. : traits<typename Transpose<typename MatrixType::PlainObject>::PlainObject>
  532. {
  533. enum { Flags = 0 };
  534. };
  535. template<typename DstXprType, typename MatrixType>
  536. struct Assignment<DstXprType, Inverse<CompleteOrthogonalDecomposition<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename CompleteOrthogonalDecomposition<MatrixType>::Scalar>, Dense2Dense>
  537. {
  538. typedef CompleteOrthogonalDecomposition<MatrixType> CodType;
  539. typedef Inverse<CodType> SrcXprType;
  540. static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename CodType::Scalar> &)
  541. {
  542. typedef Matrix<typename CodType::Scalar, CodType::RowsAtCompileTime, CodType::RowsAtCompileTime, 0, CodType::MaxRowsAtCompileTime, CodType::MaxRowsAtCompileTime> IdentityMatrixType;
  543. dst = src.nestedExpression().solve(IdentityMatrixType::Identity(src.cols(), src.cols()));
  544. }
  545. };
  546. } // end namespace internal
  547. /** \returns the matrix Q as a sequence of householder transformations */
  548. template <typename MatrixType>
  549. typename CompleteOrthogonalDecomposition<MatrixType>::HouseholderSequenceType
  550. CompleteOrthogonalDecomposition<MatrixType>::householderQ() const {
  551. return m_cpqr.householderQ();
  552. }
  553. /** \return the complete orthogonal decomposition of \c *this.
  554. *
  555. * \sa class CompleteOrthogonalDecomposition
  556. */
  557. template <typename Derived>
  558. const CompleteOrthogonalDecomposition<typename MatrixBase<Derived>::PlainObject>
  559. MatrixBase<Derived>::completeOrthogonalDecomposition() const {
  560. return CompleteOrthogonalDecomposition<PlainObject>(eval());
  561. }
  562. } // end namespace Eigen
  563. #endif // EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H