FullPivLU.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_LU_H
00011 #define EIGEN_LU_H
00012 
00013 namespace Eigen { 
00014 
00045 template<typename _MatrixType> class FullPivLU
00046 {
00047   public:
00048     typedef _MatrixType MatrixType;
00049     enum {
00050       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00051       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00052       Options = MatrixType::Options,
00053       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00054       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00055     };
00056     typedef typename MatrixType::Scalar Scalar;
00057     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00058     typedef typename internal::traits<MatrixType>::StorageKind StorageKind;
00059     typedef typename MatrixType::Index Index;
00060     typedef typename internal::plain_row_type<MatrixType, Index>::type IntRowVectorType;
00061     typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
00062     typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
00063     typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
00064 
00071     FullPivLU();
00072 
00079     FullPivLU(Index rows, Index cols);
00080 
00086     FullPivLU(const MatrixType& matrix);
00087 
00095     FullPivLU& compute(const MatrixType& matrix);
00096 
00103     inline const MatrixType& matrixLU() const
00104     {
00105       eigen_assert(m_isInitialized && "LU is not initialized.");
00106       return m_lu;
00107     }
00108 
00116     inline Index nonzeroPivots() const
00117     {
00118       eigen_assert(m_isInitialized && "LU is not initialized.");
00119       return m_nonzero_pivots;
00120     }
00121 
00125     RealScalar maxPivot() const { return m_maxpivot; }
00126 
00131     inline const PermutationPType& permutationP() const
00132     {
00133       eigen_assert(m_isInitialized && "LU is not initialized.");
00134       return m_p;
00135     }
00136 
00141     inline const PermutationQType& permutationQ() const
00142     {
00143       eigen_assert(m_isInitialized && "LU is not initialized.");
00144       return m_q;
00145     }
00146 
00161     inline const internal::kernel_retval<FullPivLU> kernel() const
00162     {
00163       eigen_assert(m_isInitialized && "LU is not initialized.");
00164       return internal::kernel_retval<FullPivLU>(*this);
00165     }
00166 
00186     inline const internal::image_retval<FullPivLU>
00187       image(const MatrixType& originalMatrix) const
00188     {
00189       eigen_assert(m_isInitialized && "LU is not initialized.");
00190       return internal::image_retval<FullPivLU>(*this, originalMatrix);
00191     }
00192 
00212     template<typename Rhs>
00213     inline const internal::solve_retval<FullPivLU, Rhs>
00214     solve(const MatrixBase<Rhs>& b) const
00215     {
00216       eigen_assert(m_isInitialized && "LU is not initialized.");
00217       return internal::solve_retval<FullPivLU, Rhs>(*this, b.derived());
00218     }
00219 
00235     typename internal::traits<MatrixType>::Scalar determinant() const;
00236 
00254     FullPivLU& setThreshold(const RealScalar& threshold)
00255     {
00256       m_usePrescribedThreshold = true;
00257       m_prescribedThreshold = threshold;
00258       return *this;
00259     }
00260 
00269     FullPivLU& setThreshold(Default_t)
00270     {
00271       m_usePrescribedThreshold = false;
00272       return *this;
00273     }
00274 
00279     RealScalar threshold() const
00280     {
00281       eigen_assert(m_isInitialized || m_usePrescribedThreshold);
00282       return m_usePrescribedThreshold ? m_prescribedThreshold
00283       // this formula comes from experimenting (see "LU precision tuning" thread on the list)
00284       // and turns out to be identical to Higham's formula used already in LDLt.
00285                                       : NumTraits<Scalar>::epsilon() * m_lu.diagonalSize();
00286     }
00287 
00294     inline Index rank() const
00295     {
00296       eigen_assert(m_isInitialized && "LU is not initialized.");
00297       RealScalar premultiplied_threshold = internal::abs(m_maxpivot) * threshold();
00298       Index result = 0;
00299       for(Index i = 0; i < m_nonzero_pivots; ++i)
00300         result += (internal::abs(m_lu.coeff(i,i)) > premultiplied_threshold);
00301       return result;
00302     }
00303 
00310     inline Index dimensionOfKernel() const
00311     {
00312       eigen_assert(m_isInitialized && "LU is not initialized.");
00313       return cols() - rank();
00314     }
00315 
00323     inline bool isInjective() const
00324     {
00325       eigen_assert(m_isInitialized && "LU is not initialized.");
00326       return rank() == cols();
00327     }
00328 
00336     inline bool isSurjective() const
00337     {
00338       eigen_assert(m_isInitialized && "LU is not initialized.");
00339       return rank() == rows();
00340     }
00341 
00348     inline bool isInvertible() const
00349     {
00350       eigen_assert(m_isInitialized && "LU is not initialized.");
00351       return isInjective() && (m_lu.rows() == m_lu.cols());
00352     }
00353 
00361     inline const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> inverse() const
00362     {
00363       eigen_assert(m_isInitialized && "LU is not initialized.");
00364       eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
00365       return internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
00366                (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
00367     }
00368 
00369     MatrixType reconstructedMatrix() const;
00370 
00371     inline Index rows() const { return m_lu.rows(); }
00372     inline Index cols() const { return m_lu.cols(); }
00373 
00374   protected:
00375     MatrixType m_lu;
00376     PermutationPType m_p;
00377     PermutationQType m_q;
00378     IntColVectorType m_rowsTranspositions;
00379     IntRowVectorType m_colsTranspositions;
00380     Index m_det_pq, m_nonzero_pivots;
00381     RealScalar m_maxpivot, m_prescribedThreshold;
00382     bool m_isInitialized, m_usePrescribedThreshold;
00383 };
00384 
00385 template<typename MatrixType>
00386 FullPivLU<MatrixType>::FullPivLU()
00387   : m_isInitialized(false), m_usePrescribedThreshold(false)
00388 {
00389 }
00390 
00391 template<typename MatrixType>
00392 FullPivLU<MatrixType>::FullPivLU(Index rows, Index cols)
00393   : m_lu(rows, cols),
00394     m_p(rows),
00395     m_q(cols),
00396     m_rowsTranspositions(rows),
00397     m_colsTranspositions(cols),
00398     m_isInitialized(false),
00399     m_usePrescribedThreshold(false)
00400 {
00401 }
00402 
00403 template<typename MatrixType>
00404 FullPivLU<MatrixType>::FullPivLU(const MatrixType& matrix)
00405   : m_lu(matrix.rows(), matrix.cols()),
00406     m_p(matrix.rows()),
00407     m_q(matrix.cols()),
00408     m_rowsTranspositions(matrix.rows()),
00409     m_colsTranspositions(matrix.cols()),
00410     m_isInitialized(false),
00411     m_usePrescribedThreshold(false)
00412 {
00413   compute(matrix);
00414 }
00415 
00416 template<typename MatrixType>
00417 FullPivLU<MatrixType>& FullPivLU<MatrixType>::compute(const MatrixType& matrix)
00418 {
00419   m_isInitialized = true;
00420   m_lu = matrix;
00421 
00422   const Index size = matrix.diagonalSize();
00423   const Index rows = matrix.rows();
00424   const Index cols = matrix.cols();
00425 
00426   // will store the transpositions, before we accumulate them at the end.
00427   // can't accumulate on-the-fly because that will be done in reverse order for the rows.
00428   m_rowsTranspositions.resize(matrix.rows());
00429   m_colsTranspositions.resize(matrix.cols());
00430   Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
00431 
00432   m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
00433   m_maxpivot = RealScalar(0);
00434 
00435   for(Index k = 0; k < size; ++k)
00436   {
00437     // First, we need to find the pivot.
00438 
00439     // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
00440     Index row_of_biggest_in_corner, col_of_biggest_in_corner;
00441     RealScalar biggest_in_corner;
00442     biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
00443                         .cwiseAbs()
00444                         .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
00445     row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
00446     col_of_biggest_in_corner += k; // need to add k to them.
00447 
00448     if(biggest_in_corner==RealScalar(0))
00449     {
00450       // before exiting, make sure to initialize the still uninitialized transpositions
00451       // in a sane state without destroying what we already have.
00452       m_nonzero_pivots = k;
00453       for(Index i = k; i < size; ++i)
00454       {
00455         m_rowsTranspositions.coeffRef(i) = i;
00456         m_colsTranspositions.coeffRef(i) = i;
00457       }
00458       break;
00459     }
00460 
00461     if(biggest_in_corner > m_maxpivot) m_maxpivot = biggest_in_corner;
00462 
00463     // Now that we've found the pivot, we need to apply the row/col swaps to
00464     // bring it to the location (k,k).
00465 
00466     m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
00467     m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
00468     if(k != row_of_biggest_in_corner) {
00469       m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
00470       ++number_of_transpositions;
00471     }
00472     if(k != col_of_biggest_in_corner) {
00473       m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
00474       ++number_of_transpositions;
00475     }
00476 
00477     // Now that the pivot is at the right location, we update the remaining
00478     // bottom-right corner by Gaussian elimination.
00479 
00480     if(k<rows-1)
00481       m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
00482     if(k<size-1)
00483       m_lu.block(k+1,k+1,rows-k-1,cols-k-1).noalias() -= m_lu.col(k).tail(rows-k-1) * m_lu.row(k).tail(cols-k-1);
00484   }
00485 
00486   // the main loop is over, we still have to accumulate the transpositions to find the
00487   // permutations P and Q
00488 
00489   m_p.setIdentity(rows);
00490   for(Index k = size-1; k >= 0; --k)
00491     m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
00492 
00493   m_q.setIdentity(cols);
00494   for(Index k = 0; k < size; ++k)
00495     m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
00496 
00497   m_det_pq = (number_of_transpositions%2) ? -1 : 1;
00498   return *this;
00499 }
00500 
00501 template<typename MatrixType>
00502 typename internal::traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() const
00503 {
00504   eigen_assert(m_isInitialized && "LU is not initialized.");
00505   eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
00506   return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
00507 }
00508 
00512 template<typename MatrixType>
00513 MatrixType FullPivLU<MatrixType>::reconstructedMatrix() const
00514 {
00515   eigen_assert(m_isInitialized && "LU is not initialized.");
00516   const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
00517   // LU
00518   MatrixType res(m_lu.rows(),m_lu.cols());
00519   // FIXME the .toDenseMatrix() should not be needed...
00520   res = m_lu.leftCols(smalldim)
00521             .template triangularView<UnitLower>().toDenseMatrix()
00522       * m_lu.topRows(smalldim)
00523             .template triangularView<Upper>().toDenseMatrix();
00524 
00525   // P^{-1}(LU)
00526   res = m_p.inverse() * res;
00527 
00528   // (P^{-1}LU)Q^{-1}
00529   res = res * m_q.inverse();
00530 
00531   return res;
00532 }
00533 
00534 /********* Implementation of kernel() **************************************************/
00535 
00536 namespace internal {
00537 template<typename _MatrixType>
00538 struct kernel_retval<FullPivLU<_MatrixType> >
00539   : kernel_retval_base<FullPivLU<_MatrixType> >
00540 {
00541   EIGEN_MAKE_KERNEL_HELPERS(FullPivLU<_MatrixType>)
00542 
00543   enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
00544             MatrixType::MaxColsAtCompileTime,
00545             MatrixType::MaxRowsAtCompileTime)
00546   };
00547 
00548   template<typename Dest> void evalTo(Dest& dst) const
00549   {
00550     const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
00551     if(dimker == 0)
00552     {
00553       // The Kernel is just {0}, so it doesn't have a basis properly speaking, but let's
00554       // avoid crashing/asserting as that depends on floating point calculations. Let's
00555       // just return a single column vector filled with zeros.
00556       dst.setZero();
00557       return;
00558     }
00559 
00560     /* Let us use the following lemma:
00561       *
00562       * Lemma: If the matrix A has the LU decomposition PAQ = LU,
00563       * then Ker A = Q(Ker U).
00564       *
00565       * Proof: trivial: just keep in mind that P, Q, L are invertible.
00566       */
00567 
00568     /* Thus, all we need to do is to compute Ker U, and then apply Q.
00569       *
00570       * U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
00571       * Thus, the diagonal of U ends with exactly
00572       * dimKer zero's. Let us use that to construct dimKer linearly
00573       * independent vectors in Ker U.
00574       */
00575 
00576     Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
00577     RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
00578     Index p = 0;
00579     for(Index i = 0; i < dec().nonzeroPivots(); ++i)
00580       if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
00581         pivots.coeffRef(p++) = i;
00582     eigen_internal_assert(p == rank());
00583 
00584     // we construct a temporaty trapezoid matrix m, by taking the U matrix and
00585     // permuting the rows and cols to bring the nonnegligible pivots to the top of
00586     // the main diagonal. We need that to be able to apply our triangular solvers.
00587     // FIXME when we get triangularView-for-rectangular-matrices, this can be simplified
00588     Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
00589            MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
00590       m(dec().matrixLU().block(0, 0, rank(), cols));
00591     for(Index i = 0; i < rank(); ++i)
00592     {
00593       if(i) m.row(i).head(i).setZero();
00594       m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
00595     }
00596     m.block(0, 0, rank(), rank());
00597     m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
00598     for(Index i = 0; i < rank(); ++i)
00599       m.col(i).swap(m.col(pivots.coeff(i)));
00600 
00601     // ok, we have our trapezoid matrix, we can apply the triangular solver.
00602     // notice that the math behind this suggests that we should apply this to the
00603     // negative of the RHS, but for performance we just put the negative sign elsewhere, see below.
00604     m.topLeftCorner(rank(), rank())
00605      .template triangularView<Upper>().solveInPlace(
00606         m.topRightCorner(rank(), dimker)
00607       );
00608 
00609     // now we must undo the column permutation that we had applied!
00610     for(Index i = rank()-1; i >= 0; --i)
00611       m.col(i).swap(m.col(pivots.coeff(i)));
00612 
00613     // see the negative sign in the next line, that's what we were talking about above.
00614     for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
00615     for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
00616     for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
00617   }
00618 };
00619 
00620 /***** Implementation of image() *****************************************************/
00621 
00622 template<typename _MatrixType>
00623 struct image_retval<FullPivLU<_MatrixType> >
00624   : image_retval_base<FullPivLU<_MatrixType> >
00625 {
00626   EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
00627 
00628   enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
00629             MatrixType::MaxColsAtCompileTime,
00630             MatrixType::MaxRowsAtCompileTime)
00631   };
00632 
00633   template<typename Dest> void evalTo(Dest& dst) const
00634   {
00635     if(rank() == 0)
00636     {
00637       // The Image is just {0}, so it doesn't have a basis properly speaking, but let's
00638       // avoid crashing/asserting as that depends on floating point calculations. Let's
00639       // just return a single column vector filled with zeros.
00640       dst.setZero();
00641       return;
00642     }
00643 
00644     Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
00645     RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
00646     Index p = 0;
00647     for(Index i = 0; i < dec().nonzeroPivots(); ++i)
00648       if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
00649         pivots.coeffRef(p++) = i;
00650     eigen_internal_assert(p == rank());
00651 
00652     for(Index i = 0; i < rank(); ++i)
00653       dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
00654   }
00655 };
00656 
00657 /***** Implementation of solve() *****************************************************/
00658 
00659 template<typename _MatrixType, typename Rhs>
00660 struct solve_retval<FullPivLU<_MatrixType>, Rhs>
00661   : solve_retval_base<FullPivLU<_MatrixType>, Rhs>
00662 {
00663   EIGEN_MAKE_SOLVE_HELPERS(FullPivLU<_MatrixType>,Rhs)
00664 
00665   template<typename Dest> void evalTo(Dest& dst) const
00666   {
00667     /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
00668      * So we proceed as follows:
00669      * Step 1: compute c = P * rhs.
00670      * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
00671      * Step 3: replace c by the solution x to Ux = c. May or may not exist.
00672      * Step 4: result = Q * c;
00673      */
00674 
00675     const Index rows = dec().rows(), cols = dec().cols(),
00676               nonzero_pivots = dec().nonzeroPivots();
00677     eigen_assert(rhs().rows() == rows);
00678     const Index smalldim = (std::min)(rows, cols);
00679 
00680     if(nonzero_pivots == 0)
00681     {
00682       dst.setZero();
00683       return;
00684     }
00685 
00686     typename Rhs::PlainObject c(rhs().rows(), rhs().cols());
00687 
00688     // Step 1
00689     c = dec().permutationP() * rhs();
00690 
00691     // Step 2
00692     dec().matrixLU()
00693         .topLeftCorner(smalldim,smalldim)
00694         .template triangularView<UnitLower>()
00695         .solveInPlace(c.topRows(smalldim));
00696     if(rows>cols)
00697     {
00698       c.bottomRows(rows-cols)
00699         -= dec().matrixLU().bottomRows(rows-cols)
00700          * c.topRows(cols);
00701     }
00702 
00703     // Step 3
00704     dec().matrixLU()
00705         .topLeftCorner(nonzero_pivots, nonzero_pivots)
00706         .template triangularView<Upper>()
00707         .solveInPlace(c.topRows(nonzero_pivots));
00708 
00709     // Step 4
00710     for(Index i = 0; i < nonzero_pivots; ++i)
00711       dst.row(dec().permutationQ().indices().coeff(i)) = c.row(i);
00712     for(Index i = nonzero_pivots; i < dec().matrixLU().cols(); ++i)
00713       dst.row(dec().permutationQ().indices().coeff(i)).setZero();
00714   }
00715 };
00716 
00717 } // end namespace internal
00718 
00719 /******* MatrixBase methods *****************************************************************/
00720 
00727 template<typename Derived>
00728 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
00729 MatrixBase<Derived>::fullPivLu() const
00730 {
00731   return FullPivLU<PlainObject>(eval());
00732 }
00733 
00734 } // end namespace Eigen
00735 
00736 #endif // EIGEN_LU_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:10:37