PartialPivLU.h 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
  5. // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
  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. #ifndef EIGEN_PARTIALLU_H
  11. #define EIGEN_PARTIALLU_H
  12. namespace Eigen {
  13. namespace internal {
  14. template<typename _MatrixType> struct traits<PartialPivLU<_MatrixType> >
  15. : traits<_MatrixType>
  16. {
  17. typedef MatrixXpr XprKind;
  18. typedef SolverStorage StorageKind;
  19. typedef int StorageIndex;
  20. typedef traits<_MatrixType> BaseTraits;
  21. enum {
  22. Flags = BaseTraits::Flags & RowMajorBit,
  23. CoeffReadCost = Dynamic
  24. };
  25. };
  26. template<typename T,typename Derived>
  27. struct enable_if_ref;
  28. // {
  29. // typedef Derived type;
  30. // };
  31. template<typename T,typename Derived>
  32. struct enable_if_ref<Ref<T>,Derived> {
  33. typedef Derived type;
  34. };
  35. } // end namespace internal
  36. /** \ingroup LU_Module
  37. *
  38. * \class PartialPivLU
  39. *
  40. * \brief LU decomposition of a matrix with partial pivoting, and related features
  41. *
  42. * \tparam _MatrixType the type of the matrix of which we are computing the LU decomposition
  43. *
  44. * This class represents a LU decomposition of a \b square \b invertible matrix, with partial pivoting: the matrix A
  45. * is decomposed as A = PLU where L is unit-lower-triangular, U is upper-triangular, and P
  46. * is a permutation matrix.
  47. *
  48. * Typically, partial pivoting LU decomposition is only considered numerically stable for square invertible
  49. * matrices. Thus LAPACK's dgesv and dgesvx require the matrix to be square and invertible. The present class
  50. * does the same. It will assert that the matrix is square, but it won't (actually it can't) check that the
  51. * matrix is invertible: it is your task to check that you only use this decomposition on invertible matrices.
  52. *
  53. * The guaranteed safe alternative, working for all matrices, is the full pivoting LU decomposition, provided
  54. * by class FullPivLU.
  55. *
  56. * This is \b not a rank-revealing LU decomposition. Many features are intentionally absent from this class,
  57. * such as rank computation. If you need these features, use class FullPivLU.
  58. *
  59. * This LU decomposition is suitable to invert invertible matrices. It is what MatrixBase::inverse() uses
  60. * in the general case.
  61. * On the other hand, it is \b not suitable to determine whether a given matrix is invertible.
  62. *
  63. * The data of the LU decomposition can be directly accessed through the methods matrixLU(), permutationP().
  64. *
  65. * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
  66. *
  67. * \sa MatrixBase::partialPivLu(), MatrixBase::determinant(), MatrixBase::inverse(), MatrixBase::computeInverse(), class FullPivLU
  68. */
  69. template<typename _MatrixType> class PartialPivLU
  70. : public SolverBase<PartialPivLU<_MatrixType> >
  71. {
  72. public:
  73. typedef _MatrixType MatrixType;
  74. typedef SolverBase<PartialPivLU> Base;
  75. friend class SolverBase<PartialPivLU>;
  76. EIGEN_GENERIC_PUBLIC_INTERFACE(PartialPivLU)
  77. enum {
  78. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  79. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  80. };
  81. typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
  82. typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
  83. typedef typename MatrixType::PlainObject PlainObject;
  84. /**
  85. * \brief Default Constructor.
  86. *
  87. * The default constructor is useful in cases in which the user intends to
  88. * perform decompositions via PartialPivLU::compute(const MatrixType&).
  89. */
  90. PartialPivLU();
  91. /** \brief Default Constructor with memory preallocation
  92. *
  93. * Like the default constructor but with preallocation of the internal data
  94. * according to the specified problem \a size.
  95. * \sa PartialPivLU()
  96. */
  97. explicit PartialPivLU(Index size);
  98. /** Constructor.
  99. *
  100. * \param matrix the matrix of which to compute the LU decomposition.
  101. *
  102. * \warning The matrix should have full rank (e.g. if it's square, it should be invertible).
  103. * If you need to deal with non-full rank, use class FullPivLU instead.
  104. */
  105. template<typename InputType>
  106. explicit PartialPivLU(const EigenBase<InputType>& matrix);
  107. /** Constructor for \link InplaceDecomposition inplace decomposition \endlink
  108. *
  109. * \param matrix the matrix of which to compute the LU decomposition.
  110. *
  111. * \warning The matrix should have full rank (e.g. if it's square, it should be invertible).
  112. * If you need to deal with non-full rank, use class FullPivLU instead.
  113. */
  114. template<typename InputType>
  115. explicit PartialPivLU(EigenBase<InputType>& matrix);
  116. template<typename InputType>
  117. PartialPivLU& compute(const EigenBase<InputType>& matrix) {
  118. m_lu = matrix.derived();
  119. compute();
  120. return *this;
  121. }
  122. /** \returns the LU decomposition matrix: the upper-triangular part is U, the
  123. * unit-lower-triangular part is L (at least for square matrices; in the non-square
  124. * case, special care is needed, see the documentation of class FullPivLU).
  125. *
  126. * \sa matrixL(), matrixU()
  127. */
  128. inline const MatrixType& matrixLU() const
  129. {
  130. eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
  131. return m_lu;
  132. }
  133. /** \returns the permutation matrix P.
  134. */
  135. inline const PermutationType& permutationP() const
  136. {
  137. eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
  138. return m_p;
  139. }
  140. #ifdef EIGEN_PARSED_BY_DOXYGEN
  141. /** This method returns the solution x to the equation Ax=b, where A is the matrix of which
  142. * *this is the LU decomposition.
  143. *
  144. * \param b the right-hand-side of the equation to solve. Can be a vector or a matrix,
  145. * the only requirement in order for the equation to make sense is that
  146. * b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.
  147. *
  148. * \returns the solution.
  149. *
  150. * Example: \include PartialPivLU_solve.cpp
  151. * Output: \verbinclude PartialPivLU_solve.out
  152. *
  153. * Since this PartialPivLU class assumes anyway that the matrix A is invertible, the solution
  154. * theoretically exists and is unique regardless of b.
  155. *
  156. * \sa TriangularView::solve(), inverse(), computeInverse()
  157. */
  158. template<typename Rhs>
  159. inline const Solve<PartialPivLU, Rhs>
  160. solve(const MatrixBase<Rhs>& b) const;
  161. #endif
  162. /** \returns an estimate of the reciprocal condition number of the matrix of which \c *this is
  163. the LU decomposition.
  164. */
  165. inline RealScalar rcond() const
  166. {
  167. eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
  168. return internal::rcond_estimate_helper(m_l1_norm, *this);
  169. }
  170. /** \returns the inverse of the matrix of which *this is the LU decomposition.
  171. *
  172. * \warning The matrix being decomposed here is assumed to be invertible. If you need to check for
  173. * invertibility, use class FullPivLU instead.
  174. *
  175. * \sa MatrixBase::inverse(), LU::inverse()
  176. */
  177. inline const Inverse<PartialPivLU> inverse() const
  178. {
  179. eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
  180. return Inverse<PartialPivLU>(*this);
  181. }
  182. /** \returns the determinant of the matrix of which
  183. * *this is the LU decomposition. It has only linear complexity
  184. * (that is, O(n) where n is the dimension of the square matrix)
  185. * as the LU decomposition has already been computed.
  186. *
  187. * \note For fixed-size matrices of size up to 4, MatrixBase::determinant() offers
  188. * optimized paths.
  189. *
  190. * \warning a determinant can be very big or small, so for matrices
  191. * of large enough dimension, there is a risk of overflow/underflow.
  192. *
  193. * \sa MatrixBase::determinant()
  194. */
  195. Scalar determinant() const;
  196. MatrixType reconstructedMatrix() const;
  197. EIGEN_CONSTEXPR inline Index rows() const EIGEN_NOEXCEPT { return m_lu.rows(); }
  198. EIGEN_CONSTEXPR inline Index cols() const EIGEN_NOEXCEPT { return m_lu.cols(); }
  199. #ifndef EIGEN_PARSED_BY_DOXYGEN
  200. template<typename RhsType, typename DstType>
  201. EIGEN_DEVICE_FUNC
  202. void _solve_impl(const RhsType &rhs, DstType &dst) const {
  203. /* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
  204. * So we proceed as follows:
  205. * Step 1: compute c = Pb.
  206. * Step 2: replace c by the solution x to Lx = c.
  207. * Step 3: replace c by the solution x to Ux = c.
  208. */
  209. // Step 1
  210. dst = permutationP() * rhs;
  211. // Step 2
  212. m_lu.template triangularView<UnitLower>().solveInPlace(dst);
  213. // Step 3
  214. m_lu.template triangularView<Upper>().solveInPlace(dst);
  215. }
  216. template<bool Conjugate, typename RhsType, typename DstType>
  217. EIGEN_DEVICE_FUNC
  218. void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const {
  219. /* The decomposition PA = LU can be rewritten as A^T = U^T L^T P.
  220. * So we proceed as follows:
  221. * Step 1: compute c as the solution to L^T c = b
  222. * Step 2: replace c by the solution x to U^T x = c.
  223. * Step 3: update c = P^-1 c.
  224. */
  225. eigen_assert(rhs.rows() == m_lu.cols());
  226. // Step 1
  227. dst = m_lu.template triangularView<Upper>().transpose()
  228. .template conjugateIf<Conjugate>().solve(rhs);
  229. // Step 2
  230. m_lu.template triangularView<UnitLower>().transpose()
  231. .template conjugateIf<Conjugate>().solveInPlace(dst);
  232. // Step 3
  233. dst = permutationP().transpose() * dst;
  234. }
  235. #endif
  236. protected:
  237. static void check_template_parameters()
  238. {
  239. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  240. }
  241. void compute();
  242. MatrixType m_lu;
  243. PermutationType m_p;
  244. TranspositionType m_rowsTranspositions;
  245. RealScalar m_l1_norm;
  246. signed char m_det_p;
  247. bool m_isInitialized;
  248. };
  249. template<typename MatrixType>
  250. PartialPivLU<MatrixType>::PartialPivLU()
  251. : m_lu(),
  252. m_p(),
  253. m_rowsTranspositions(),
  254. m_l1_norm(0),
  255. m_det_p(0),
  256. m_isInitialized(false)
  257. {
  258. }
  259. template<typename MatrixType>
  260. PartialPivLU<MatrixType>::PartialPivLU(Index size)
  261. : m_lu(size, size),
  262. m_p(size),
  263. m_rowsTranspositions(size),
  264. m_l1_norm(0),
  265. m_det_p(0),
  266. m_isInitialized(false)
  267. {
  268. }
  269. template<typename MatrixType>
  270. template<typename InputType>
  271. PartialPivLU<MatrixType>::PartialPivLU(const EigenBase<InputType>& matrix)
  272. : m_lu(matrix.rows(),matrix.cols()),
  273. m_p(matrix.rows()),
  274. m_rowsTranspositions(matrix.rows()),
  275. m_l1_norm(0),
  276. m_det_p(0),
  277. m_isInitialized(false)
  278. {
  279. compute(matrix.derived());
  280. }
  281. template<typename MatrixType>
  282. template<typename InputType>
  283. PartialPivLU<MatrixType>::PartialPivLU(EigenBase<InputType>& matrix)
  284. : m_lu(matrix.derived()),
  285. m_p(matrix.rows()),
  286. m_rowsTranspositions(matrix.rows()),
  287. m_l1_norm(0),
  288. m_det_p(0),
  289. m_isInitialized(false)
  290. {
  291. compute();
  292. }
  293. namespace internal {
  294. /** \internal This is the blocked version of fullpivlu_unblocked() */
  295. template<typename Scalar, int StorageOrder, typename PivIndex, int SizeAtCompileTime=Dynamic>
  296. struct partial_lu_impl
  297. {
  298. static const int UnBlockedBound = 16;
  299. static const bool UnBlockedAtCompileTime = SizeAtCompileTime!=Dynamic && SizeAtCompileTime<=UnBlockedBound;
  300. static const int ActualSizeAtCompileTime = UnBlockedAtCompileTime ? SizeAtCompileTime : Dynamic;
  301. // Remaining rows and columns at compile-time:
  302. static const int RRows = SizeAtCompileTime==2 ? 1 : Dynamic;
  303. static const int RCols = SizeAtCompileTime==2 ? 1 : Dynamic;
  304. typedef Matrix<Scalar, ActualSizeAtCompileTime, ActualSizeAtCompileTime, StorageOrder> MatrixType;
  305. typedef Ref<MatrixType> MatrixTypeRef;
  306. typedef Ref<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > BlockType;
  307. typedef typename MatrixType::RealScalar RealScalar;
  308. /** \internal performs the LU decomposition in-place of the matrix \a lu
  309. * using an unblocked algorithm.
  310. *
  311. * In addition, this function returns the row transpositions in the
  312. * vector \a row_transpositions which must have a size equal to the number
  313. * of columns of the matrix \a lu, and an integer \a nb_transpositions
  314. * which returns the actual number of transpositions.
  315. *
  316. * \returns The index of the first pivot which is exactly zero if any, or a negative number otherwise.
  317. */
  318. static Index unblocked_lu(MatrixTypeRef& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
  319. {
  320. typedef scalar_score_coeff_op<Scalar> Scoring;
  321. typedef typename Scoring::result_type Score;
  322. const Index rows = lu.rows();
  323. const Index cols = lu.cols();
  324. const Index size = (std::min)(rows,cols);
  325. // For small compile-time matrices it is worth processing the last row separately:
  326. // speedup: +100% for 2x2, +10% for others.
  327. const Index endk = UnBlockedAtCompileTime ? size-1 : size;
  328. nb_transpositions = 0;
  329. Index first_zero_pivot = -1;
  330. for(Index k = 0; k < endk; ++k)
  331. {
  332. int rrows = internal::convert_index<int>(rows-k-1);
  333. int rcols = internal::convert_index<int>(cols-k-1);
  334. Index row_of_biggest_in_col;
  335. Score biggest_in_corner
  336. = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
  337. row_of_biggest_in_col += k;
  338. row_transpositions[k] = PivIndex(row_of_biggest_in_col);
  339. if(biggest_in_corner != Score(0))
  340. {
  341. if(k != row_of_biggest_in_col)
  342. {
  343. lu.row(k).swap(lu.row(row_of_biggest_in_col));
  344. ++nb_transpositions;
  345. }
  346. lu.col(k).tail(fix<RRows>(rrows)) /= lu.coeff(k,k);
  347. }
  348. else if(first_zero_pivot==-1)
  349. {
  350. // the pivot is exactly zero, we record the index of the first pivot which is exactly 0,
  351. // and continue the factorization such we still have A = PLU
  352. first_zero_pivot = k;
  353. }
  354. if(k<rows-1)
  355. lu.bottomRightCorner(fix<RRows>(rrows),fix<RCols>(rcols)).noalias() -= lu.col(k).tail(fix<RRows>(rrows)) * lu.row(k).tail(fix<RCols>(rcols));
  356. }
  357. // special handling of the last entry
  358. if(UnBlockedAtCompileTime)
  359. {
  360. Index k = endk;
  361. row_transpositions[k] = PivIndex(k);
  362. if (Scoring()(lu(k, k)) == Score(0) && first_zero_pivot == -1)
  363. first_zero_pivot = k;
  364. }
  365. return first_zero_pivot;
  366. }
  367. /** \internal performs the LU decomposition in-place of the matrix represented
  368. * by the variables \a rows, \a cols, \a lu_data, and \a lu_stride using a
  369. * recursive, blocked algorithm.
  370. *
  371. * In addition, this function returns the row transpositions in the
  372. * vector \a row_transpositions which must have a size equal to the number
  373. * of columns of the matrix \a lu, and an integer \a nb_transpositions
  374. * which returns the actual number of transpositions.
  375. *
  376. * \returns The index of the first pivot which is exactly zero if any, or a negative number otherwise.
  377. *
  378. * \note This very low level interface using pointers, etc. is to:
  379. * 1 - reduce the number of instantiations to the strict minimum
  380. * 2 - avoid infinite recursion of the instantiations with Block<Block<Block<...> > >
  381. */
  382. static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256)
  383. {
  384. MatrixTypeRef lu = MatrixType::Map(lu_data,rows, cols, OuterStride<>(luStride));
  385. const Index size = (std::min)(rows,cols);
  386. // if the matrix is too small, no blocking:
  387. if(UnBlockedAtCompileTime || size<=UnBlockedBound)
  388. {
  389. return unblocked_lu(lu, row_transpositions, nb_transpositions);
  390. }
  391. // automatically adjust the number of subdivisions to the size
  392. // of the matrix so that there is enough sub blocks:
  393. Index blockSize;
  394. {
  395. blockSize = size/8;
  396. blockSize = (blockSize/16)*16;
  397. blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize);
  398. }
  399. nb_transpositions = 0;
  400. Index first_zero_pivot = -1;
  401. for(Index k = 0; k < size; k+=blockSize)
  402. {
  403. Index bs = (std::min)(size-k,blockSize); // actual size of the block
  404. Index trows = rows - k - bs; // trailing rows
  405. Index tsize = size - k - bs; // trailing size
  406. // partition the matrix:
  407. // A00 | A01 | A02
  408. // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
  409. // A20 | A21 | A22
  410. BlockType A_0 = lu.block(0,0,rows,k);
  411. BlockType A_2 = lu.block(0,k+bs,rows,tsize);
  412. BlockType A11 = lu.block(k,k,bs,bs);
  413. BlockType A12 = lu.block(k,k+bs,bs,tsize);
  414. BlockType A21 = lu.block(k+bs,k,trows,bs);
  415. BlockType A22 = lu.block(k+bs,k+bs,trows,tsize);
  416. PivIndex nb_transpositions_in_panel;
  417. // recursively call the blocked LU algorithm on [A11^T A21^T]^T
  418. // with a very small blocking size:
  419. Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
  420. row_transpositions+k, nb_transpositions_in_panel, 16);
  421. if(ret>=0 && first_zero_pivot==-1)
  422. first_zero_pivot = k+ret;
  423. nb_transpositions += nb_transpositions_in_panel;
  424. // update permutations and apply them to A_0
  425. for(Index i=k; i<k+bs; ++i)
  426. {
  427. Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
  428. A_0.row(i).swap(A_0.row(piv));
  429. }
  430. if(trows)
  431. {
  432. // apply permutations to A_2
  433. for(Index i=k;i<k+bs; ++i)
  434. A_2.row(i).swap(A_2.row(row_transpositions[i]));
  435. // A12 = A11^-1 A12
  436. A11.template triangularView<UnitLower>().solveInPlace(A12);
  437. A22.noalias() -= A21 * A12;
  438. }
  439. }
  440. return first_zero_pivot;
  441. }
  442. };
  443. /** \internal performs the LU decomposition with partial pivoting in-place.
  444. */
  445. template<typename MatrixType, typename TranspositionType>
  446. void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename TranspositionType::StorageIndex& nb_transpositions)
  447. {
  448. // Special-case of zero matrix.
  449. if (lu.rows() == 0 || lu.cols() == 0) {
  450. nb_transpositions = 0;
  451. return;
  452. }
  453. eigen_assert(lu.cols() == row_transpositions.size());
  454. eigen_assert(row_transpositions.size() < 2 || (&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
  455. partial_lu_impl
  456. < typename MatrixType::Scalar, MatrixType::Flags&RowMajorBit?RowMajor:ColMajor,
  457. typename TranspositionType::StorageIndex,
  458. EIGEN_SIZE_MIN_PREFER_FIXED(MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime)>
  459. ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
  460. }
  461. } // end namespace internal
  462. template<typename MatrixType>
  463. void PartialPivLU<MatrixType>::compute()
  464. {
  465. check_template_parameters();
  466. // the row permutation is stored as int indices, so just to be sure:
  467. eigen_assert(m_lu.rows()<NumTraits<int>::highest());
  468. if(m_lu.cols()>0)
  469. m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
  470. else
  471. m_l1_norm = RealScalar(0);
  472. eigen_assert(m_lu.rows() == m_lu.cols() && "PartialPivLU is only for square (and moreover invertible) matrices");
  473. const Index size = m_lu.rows();
  474. m_rowsTranspositions.resize(size);
  475. typename TranspositionType::StorageIndex nb_transpositions;
  476. internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
  477. m_det_p = (nb_transpositions%2) ? -1 : 1;
  478. m_p = m_rowsTranspositions;
  479. m_isInitialized = true;
  480. }
  481. template<typename MatrixType>
  482. typename PartialPivLU<MatrixType>::Scalar PartialPivLU<MatrixType>::determinant() const
  483. {
  484. eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
  485. return Scalar(m_det_p) * m_lu.diagonal().prod();
  486. }
  487. /** \returns the matrix represented by the decomposition,
  488. * i.e., it returns the product: P^{-1} L U.
  489. * This function is provided for debug purpose. */
  490. template<typename MatrixType>
  491. MatrixType PartialPivLU<MatrixType>::reconstructedMatrix() const
  492. {
  493. eigen_assert(m_isInitialized && "LU is not initialized.");
  494. // LU
  495. MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
  496. * m_lu.template triangularView<Upper>();
  497. // P^{-1}(LU)
  498. res = m_p.inverse() * res;
  499. return res;
  500. }
  501. /***** Implementation details *****************************************************/
  502. namespace internal {
  503. /***** Implementation of inverse() *****************************************************/
  504. template<typename DstXprType, typename MatrixType>
  505. struct Assignment<DstXprType, Inverse<PartialPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename PartialPivLU<MatrixType>::Scalar>, Dense2Dense>
  506. {
  507. typedef PartialPivLU<MatrixType> LuType;
  508. typedef Inverse<LuType> SrcXprType;
  509. static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename LuType::Scalar> &)
  510. {
  511. dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
  512. }
  513. };
  514. } // end namespace internal
  515. /******** MatrixBase methods *******/
  516. /** \lu_module
  517. *
  518. * \return the partial-pivoting LU decomposition of \c *this.
  519. *
  520. * \sa class PartialPivLU
  521. */
  522. template<typename Derived>
  523. inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
  524. MatrixBase<Derived>::partialPivLu() const
  525. {
  526. return PartialPivLU<PlainObject>(eval());
  527. }
  528. /** \lu_module
  529. *
  530. * Synonym of partialPivLu().
  531. *
  532. * \return the partial-pivoting LU decomposition of \c *this.
  533. *
  534. * \sa class PartialPivLU
  535. */
  536. template<typename Derived>
  537. inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
  538. MatrixBase<Derived>::lu() const
  539. {
  540. return PartialPivLU<PlainObject>(eval());
  541. }
  542. } // end namespace Eigen
  543. #endif // EIGEN_PARTIALLU_H