Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/LU/FullPivLU.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_LU_H
0011 #define EIGEN_LU_H
0012 
0013 namespace Eigen {
0014 
0015 namespace internal {
0016 template<typename _MatrixType> struct traits<FullPivLU<_MatrixType> >
0017  : traits<_MatrixType>
0018 {
0019   typedef MatrixXpr XprKind;
0020   typedef SolverStorage StorageKind;
0021   typedef int StorageIndex;
0022   enum { Flags = 0 };
0023 };
0024 
0025 } // end namespace internal
0026 
0027 /** \ingroup LU_Module
0028   *
0029   * \class FullPivLU
0030   *
0031   * \brief LU decomposition of a matrix with complete pivoting, and related features
0032   *
0033   * \tparam _MatrixType the type of the matrix of which we are computing the LU decomposition
0034   *
0035   * This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A is
0036   * decomposed as \f$ A = P^{-1} L U Q^{-1} \f$ where L is unit-lower-triangular, U is
0037   * upper-triangular, and P and Q are permutation matrices. This is a rank-revealing LU
0038   * decomposition. The eigenvalues (diagonal coefficients) of U are sorted in such a way that any
0039   * zeros are at the end.
0040   *
0041   * This decomposition provides the generic approach to solving systems of linear equations, computing
0042   * the rank, invertibility, inverse, kernel, and determinant.
0043   *
0044   * This LU decomposition is very stable and well tested with large matrices. However there are use cases where the SVD
0045   * decomposition is inherently more stable and/or flexible. For example, when computing the kernel of a matrix,
0046   * working with the SVD allows to select the smallest singular values of the matrix, something that
0047   * the LU decomposition doesn't see.
0048   *
0049   * The data of the LU decomposition can be directly accessed through the methods matrixLU(),
0050   * permutationP(), permutationQ().
0051   *
0052   * As an example, here is how the original matrix can be retrieved:
0053   * \include class_FullPivLU.cpp
0054   * Output: \verbinclude class_FullPivLU.out
0055   *
0056   * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
0057   *
0058   * \sa MatrixBase::fullPivLu(), MatrixBase::determinant(), MatrixBase::inverse()
0059   */
0060 template<typename _MatrixType> class FullPivLU
0061   : public SolverBase<FullPivLU<_MatrixType> >
0062 {
0063   public:
0064     typedef _MatrixType MatrixType;
0065     typedef SolverBase<FullPivLU> Base;
0066     friend class SolverBase<FullPivLU>;
0067 
0068     EIGEN_GENERIC_PUBLIC_INTERFACE(FullPivLU)
0069     enum {
0070       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0071       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
0072     };
0073     typedef typename internal::plain_row_type<MatrixType, StorageIndex>::type IntRowVectorType;
0074     typedef typename internal::plain_col_type<MatrixType, StorageIndex>::type IntColVectorType;
0075     typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
0076     typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
0077     typedef typename MatrixType::PlainObject PlainObject;
0078 
0079     /**
0080       * \brief Default Constructor.
0081       *
0082       * The default constructor is useful in cases in which the user intends to
0083       * perform decompositions via LU::compute(const MatrixType&).
0084       */
0085     FullPivLU();
0086 
0087     /** \brief Default Constructor with memory preallocation
0088       *
0089       * Like the default constructor but with preallocation of the internal data
0090       * according to the specified problem \a size.
0091       * \sa FullPivLU()
0092       */
0093     FullPivLU(Index rows, Index cols);
0094 
0095     /** Constructor.
0096       *
0097       * \param matrix the matrix of which to compute the LU decomposition.
0098       *               It is required to be nonzero.
0099       */
0100     template<typename InputType>
0101     explicit FullPivLU(const EigenBase<InputType>& matrix);
0102 
0103     /** \brief Constructs a LU factorization from a given matrix
0104       *
0105       * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref.
0106       *
0107       * \sa FullPivLU(const EigenBase&)
0108       */
0109     template<typename InputType>
0110     explicit FullPivLU(EigenBase<InputType>& matrix);
0111 
0112     /** Computes the LU decomposition of the given matrix.
0113       *
0114       * \param matrix the matrix of which to compute the LU decomposition.
0115       *               It is required to be nonzero.
0116       *
0117       * \returns a reference to *this
0118       */
0119     template<typename InputType>
0120     FullPivLU& compute(const EigenBase<InputType>& matrix) {
0121       m_lu = matrix.derived();
0122       computeInPlace();
0123       return *this;
0124     }
0125 
0126     /** \returns the LU decomposition matrix: the upper-triangular part is U, the
0127       * unit-lower-triangular part is L (at least for square matrices; in the non-square
0128       * case, special care is needed, see the documentation of class FullPivLU).
0129       *
0130       * \sa matrixL(), matrixU()
0131       */
0132     inline const MatrixType& matrixLU() const
0133     {
0134       eigen_assert(m_isInitialized && "LU is not initialized.");
0135       return m_lu;
0136     }
0137 
0138     /** \returns the number of nonzero pivots in the LU decomposition.
0139       * Here nonzero is meant in the exact sense, not in a fuzzy sense.
0140       * So that notion isn't really intrinsically interesting, but it is
0141       * still useful when implementing algorithms.
0142       *
0143       * \sa rank()
0144       */
0145     inline Index nonzeroPivots() const
0146     {
0147       eigen_assert(m_isInitialized && "LU is not initialized.");
0148       return m_nonzero_pivots;
0149     }
0150 
0151     /** \returns the absolute value of the biggest pivot, i.e. the biggest
0152       *          diagonal coefficient of U.
0153       */
0154     RealScalar maxPivot() const { return m_maxpivot; }
0155 
0156     /** \returns the permutation matrix P
0157       *
0158       * \sa permutationQ()
0159       */
0160     EIGEN_DEVICE_FUNC inline const PermutationPType& permutationP() const
0161     {
0162       eigen_assert(m_isInitialized && "LU is not initialized.");
0163       return m_p;
0164     }
0165 
0166     /** \returns the permutation matrix Q
0167       *
0168       * \sa permutationP()
0169       */
0170     inline const PermutationQType& permutationQ() const
0171     {
0172       eigen_assert(m_isInitialized && "LU is not initialized.");
0173       return m_q;
0174     }
0175 
0176     /** \returns the kernel of the matrix, also called its null-space. The columns of the returned matrix
0177       * will form a basis of the kernel.
0178       *
0179       * \note If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.
0180       *
0181       * \note This method has to determine which pivots should be considered nonzero.
0182       *       For that, it uses the threshold value that you can control by calling
0183       *       setThreshold(const RealScalar&).
0184       *
0185       * Example: \include FullPivLU_kernel.cpp
0186       * Output: \verbinclude FullPivLU_kernel.out
0187       *
0188       * \sa image()
0189       */
0190     inline const internal::kernel_retval<FullPivLU> kernel() const
0191     {
0192       eigen_assert(m_isInitialized && "LU is not initialized.");
0193       return internal::kernel_retval<FullPivLU>(*this);
0194     }
0195 
0196     /** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
0197       * will form a basis of the image (column-space).
0198       *
0199       * \param originalMatrix the original matrix, of which *this is the LU decomposition.
0200       *                       The reason why it is needed to pass it here, is that this allows
0201       *                       a large optimization, as otherwise this method would need to reconstruct it
0202       *                       from the LU decomposition.
0203       *
0204       * \note If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.
0205       *
0206       * \note This method has to determine which pivots should be considered nonzero.
0207       *       For that, it uses the threshold value that you can control by calling
0208       *       setThreshold(const RealScalar&).
0209       *
0210       * Example: \include FullPivLU_image.cpp
0211       * Output: \verbinclude FullPivLU_image.out
0212       *
0213       * \sa kernel()
0214       */
0215     inline const internal::image_retval<FullPivLU>
0216       image(const MatrixType& originalMatrix) const
0217     {
0218       eigen_assert(m_isInitialized && "LU is not initialized.");
0219       return internal::image_retval<FullPivLU>(*this, originalMatrix);
0220     }
0221 
0222     #ifdef EIGEN_PARSED_BY_DOXYGEN
0223     /** \return a solution x to the equation Ax=b, where A is the matrix of which
0224       * *this is the LU decomposition.
0225       *
0226       * \param b the right-hand-side of the equation to solve. Can be a vector or a matrix,
0227       *          the only requirement in order for the equation to make sense is that
0228       *          b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.
0229       *
0230       * \returns a solution.
0231       *
0232       * \note_about_checking_solutions
0233       *
0234       * \note_about_arbitrary_choice_of_solution
0235       * \note_about_using_kernel_to_study_multiple_solutions
0236       *
0237       * Example: \include FullPivLU_solve.cpp
0238       * Output: \verbinclude FullPivLU_solve.out
0239       *
0240       * \sa TriangularView::solve(), kernel(), inverse()
0241       */
0242     template<typename Rhs>
0243     inline const Solve<FullPivLU, Rhs>
0244     solve(const MatrixBase<Rhs>& b) const;
0245     #endif
0246 
0247     /** \returns an estimate of the reciprocal condition number of the matrix of which \c *this is
0248         the LU decomposition.
0249       */
0250     inline RealScalar rcond() const
0251     {
0252       eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
0253       return internal::rcond_estimate_helper(m_l1_norm, *this);
0254     }
0255 
0256     /** \returns the determinant of the matrix of which
0257       * *this is the LU decomposition. It has only linear complexity
0258       * (that is, O(n) where n is the dimension of the square matrix)
0259       * as the LU decomposition has already been computed.
0260       *
0261       * \note This is only for square matrices.
0262       *
0263       * \note For fixed-size matrices of size up to 4, MatrixBase::determinant() offers
0264       *       optimized paths.
0265       *
0266       * \warning a determinant can be very big or small, so for matrices
0267       * of large enough dimension, there is a risk of overflow/underflow.
0268       *
0269       * \sa MatrixBase::determinant()
0270       */
0271     typename internal::traits<MatrixType>::Scalar determinant() const;
0272 
0273     /** Allows to prescribe a threshold to be used by certain methods, such as rank(),
0274       * who need to determine when pivots are to be considered nonzero. This is not used for the
0275       * LU decomposition itself.
0276       *
0277       * When it needs to get the threshold value, Eigen calls threshold(). By default, this
0278       * uses a formula to automatically determine a reasonable threshold.
0279       * Once you have called the present method setThreshold(const RealScalar&),
0280       * your value is used instead.
0281       *
0282       * \param threshold The new value to use as the threshold.
0283       *
0284       * A pivot will be considered nonzero if its absolute value is strictly greater than
0285       *  \f$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \f$
0286       * where maxpivot is the biggest pivot.
0287       *
0288       * If you want to come back to the default behavior, call setThreshold(Default_t)
0289       */
0290     FullPivLU& setThreshold(const RealScalar& threshold)
0291     {
0292       m_usePrescribedThreshold = true;
0293       m_prescribedThreshold = threshold;
0294       return *this;
0295     }
0296 
0297     /** Allows to come back to the default behavior, letting Eigen use its default formula for
0298       * determining the threshold.
0299       *
0300       * You should pass the special object Eigen::Default as parameter here.
0301       * \code lu.setThreshold(Eigen::Default); \endcode
0302       *
0303       * See the documentation of setThreshold(const RealScalar&).
0304       */
0305     FullPivLU& setThreshold(Default_t)
0306     {
0307       m_usePrescribedThreshold = false;
0308       return *this;
0309     }
0310 
0311     /** Returns the threshold that will be used by certain methods such as rank().
0312       *
0313       * See the documentation of setThreshold(const RealScalar&).
0314       */
0315     RealScalar threshold() const
0316     {
0317       eigen_assert(m_isInitialized || m_usePrescribedThreshold);
0318       return m_usePrescribedThreshold ? m_prescribedThreshold
0319       // this formula comes from experimenting (see "LU precision tuning" thread on the list)
0320       // and turns out to be identical to Higham's formula used already in LDLt.
0321           : NumTraits<Scalar>::epsilon() * RealScalar(m_lu.diagonalSize());
0322     }
0323 
0324     /** \returns the rank of the matrix of which *this is the LU decomposition.
0325       *
0326       * \note This method has to determine which pivots should be considered nonzero.
0327       *       For that, it uses the threshold value that you can control by calling
0328       *       setThreshold(const RealScalar&).
0329       */
0330     inline Index rank() const
0331     {
0332       using std::abs;
0333       eigen_assert(m_isInitialized && "LU is not initialized.");
0334       RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
0335       Index result = 0;
0336       for(Index i = 0; i < m_nonzero_pivots; ++i)
0337         result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
0338       return result;
0339     }
0340 
0341     /** \returns the dimension of the kernel of the matrix of which *this is the LU decomposition.
0342       *
0343       * \note This method has to determine which pivots should be considered nonzero.
0344       *       For that, it uses the threshold value that you can control by calling
0345       *       setThreshold(const RealScalar&).
0346       */
0347     inline Index dimensionOfKernel() const
0348     {
0349       eigen_assert(m_isInitialized && "LU is not initialized.");
0350       return cols() - rank();
0351     }
0352 
0353     /** \returns true if the matrix of which *this is the LU decomposition represents an injective
0354       *          linear map, i.e. has trivial kernel; false otherwise.
0355       *
0356       * \note This method has to determine which pivots should be considered nonzero.
0357       *       For that, it uses the threshold value that you can control by calling
0358       *       setThreshold(const RealScalar&).
0359       */
0360     inline bool isInjective() const
0361     {
0362       eigen_assert(m_isInitialized && "LU is not initialized.");
0363       return rank() == cols();
0364     }
0365 
0366     /** \returns true if the matrix of which *this is the LU decomposition represents a surjective
0367       *          linear map; false otherwise.
0368       *
0369       * \note This method has to determine which pivots should be considered nonzero.
0370       *       For that, it uses the threshold value that you can control by calling
0371       *       setThreshold(const RealScalar&).
0372       */
0373     inline bool isSurjective() const
0374     {
0375       eigen_assert(m_isInitialized && "LU is not initialized.");
0376       return rank() == rows();
0377     }
0378 
0379     /** \returns true if the matrix of which *this is the LU decomposition is invertible.
0380       *
0381       * \note This method has to determine which pivots should be considered nonzero.
0382       *       For that, it uses the threshold value that you can control by calling
0383       *       setThreshold(const RealScalar&).
0384       */
0385     inline bool isInvertible() const
0386     {
0387       eigen_assert(m_isInitialized && "LU is not initialized.");
0388       return isInjective() && (m_lu.rows() == m_lu.cols());
0389     }
0390 
0391     /** \returns the inverse of the matrix of which *this is the LU decomposition.
0392       *
0393       * \note If this matrix is not invertible, the returned matrix has undefined coefficients.
0394       *       Use isInvertible() to first determine whether this matrix is invertible.
0395       *
0396       * \sa MatrixBase::inverse()
0397       */
0398     inline const Inverse<FullPivLU> inverse() const
0399     {
0400       eigen_assert(m_isInitialized && "LU is not initialized.");
0401       eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
0402       return Inverse<FullPivLU>(*this);
0403     }
0404 
0405     MatrixType reconstructedMatrix() const;
0406 
0407     EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0408     inline Index rows() const EIGEN_NOEXCEPT { return m_lu.rows(); }
0409     EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0410     inline Index cols() const EIGEN_NOEXCEPT { return m_lu.cols(); }
0411 
0412     #ifndef EIGEN_PARSED_BY_DOXYGEN
0413     template<typename RhsType, typename DstType>
0414     void _solve_impl(const RhsType &rhs, DstType &dst) const;
0415 
0416     template<bool Conjugate, typename RhsType, typename DstType>
0417     void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
0418     #endif
0419 
0420   protected:
0421 
0422     static void check_template_parameters()
0423     {
0424       EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
0425     }
0426 
0427     void computeInPlace();
0428 
0429     MatrixType m_lu;
0430     PermutationPType m_p;
0431     PermutationQType m_q;
0432     IntColVectorType m_rowsTranspositions;
0433     IntRowVectorType m_colsTranspositions;
0434     Index m_nonzero_pivots;
0435     RealScalar m_l1_norm;
0436     RealScalar m_maxpivot, m_prescribedThreshold;
0437     signed char m_det_pq;
0438     bool m_isInitialized, m_usePrescribedThreshold;
0439 };
0440 
0441 template<typename MatrixType>
0442 FullPivLU<MatrixType>::FullPivLU()
0443   : m_isInitialized(false), m_usePrescribedThreshold(false)
0444 {
0445 }
0446 
0447 template<typename MatrixType>
0448 FullPivLU<MatrixType>::FullPivLU(Index rows, Index cols)
0449   : m_lu(rows, cols),
0450     m_p(rows),
0451     m_q(cols),
0452     m_rowsTranspositions(rows),
0453     m_colsTranspositions(cols),
0454     m_isInitialized(false),
0455     m_usePrescribedThreshold(false)
0456 {
0457 }
0458 
0459 template<typename MatrixType>
0460 template<typename InputType>
0461 FullPivLU<MatrixType>::FullPivLU(const EigenBase<InputType>& matrix)
0462   : m_lu(matrix.rows(), matrix.cols()),
0463     m_p(matrix.rows()),
0464     m_q(matrix.cols()),
0465     m_rowsTranspositions(matrix.rows()),
0466     m_colsTranspositions(matrix.cols()),
0467     m_isInitialized(false),
0468     m_usePrescribedThreshold(false)
0469 {
0470   compute(matrix.derived());
0471 }
0472 
0473 template<typename MatrixType>
0474 template<typename InputType>
0475 FullPivLU<MatrixType>::FullPivLU(EigenBase<InputType>& matrix)
0476   : m_lu(matrix.derived()),
0477     m_p(matrix.rows()),
0478     m_q(matrix.cols()),
0479     m_rowsTranspositions(matrix.rows()),
0480     m_colsTranspositions(matrix.cols()),
0481     m_isInitialized(false),
0482     m_usePrescribedThreshold(false)
0483 {
0484   computeInPlace();
0485 }
0486 
0487 template<typename MatrixType>
0488 void FullPivLU<MatrixType>::computeInPlace()
0489 {
0490   check_template_parameters();
0491 
0492   // the permutations are stored as int indices, so just to be sure:
0493   eigen_assert(m_lu.rows()<=NumTraits<int>::highest() && m_lu.cols()<=NumTraits<int>::highest());
0494 
0495   m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
0496 
0497   const Index size = m_lu.diagonalSize();
0498   const Index rows = m_lu.rows();
0499   const Index cols = m_lu.cols();
0500 
0501   // will store the transpositions, before we accumulate them at the end.
0502   // can't accumulate on-the-fly because that will be done in reverse order for the rows.
0503   m_rowsTranspositions.resize(m_lu.rows());
0504   m_colsTranspositions.resize(m_lu.cols());
0505   Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
0506 
0507   m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
0508   m_maxpivot = RealScalar(0);
0509 
0510   for(Index k = 0; k < size; ++k)
0511   {
0512     // First, we need to find the pivot.
0513 
0514     // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
0515     Index row_of_biggest_in_corner, col_of_biggest_in_corner;
0516     typedef internal::scalar_score_coeff_op<Scalar> Scoring;
0517     typedef typename Scoring::result_type Score;
0518     Score biggest_in_corner;
0519     biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
0520                         .unaryExpr(Scoring())
0521                         .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
0522     row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
0523     col_of_biggest_in_corner += k; // need to add k to them.
0524 
0525     if(biggest_in_corner==Score(0))
0526     {
0527       // before exiting, make sure to initialize the still uninitialized transpositions
0528       // in a sane state without destroying what we already have.
0529       m_nonzero_pivots = k;
0530       for(Index i = k; i < size; ++i)
0531       {
0532         m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
0533         m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
0534       }
0535       break;
0536     }
0537 
0538     RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
0539     if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
0540 
0541     // Now that we've found the pivot, we need to apply the row/col swaps to
0542     // bring it to the location (k,k).
0543 
0544     m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
0545     m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
0546     if(k != row_of_biggest_in_corner) {
0547       m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
0548       ++number_of_transpositions;
0549     }
0550     if(k != col_of_biggest_in_corner) {
0551       m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
0552       ++number_of_transpositions;
0553     }
0554 
0555     // Now that the pivot is at the right location, we update the remaining
0556     // bottom-right corner by Gaussian elimination.
0557 
0558     if(k<rows-1)
0559       m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
0560     if(k<size-1)
0561       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);
0562   }
0563 
0564   // the main loop is over, we still have to accumulate the transpositions to find the
0565   // permutations P and Q
0566 
0567   m_p.setIdentity(rows);
0568   for(Index k = size-1; k >= 0; --k)
0569     m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
0570 
0571   m_q.setIdentity(cols);
0572   for(Index k = 0; k < size; ++k)
0573     m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
0574 
0575   m_det_pq = (number_of_transpositions%2) ? -1 : 1;
0576 
0577   m_isInitialized = true;
0578 }
0579 
0580 template<typename MatrixType>
0581 typename internal::traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() const
0582 {
0583   eigen_assert(m_isInitialized && "LU is not initialized.");
0584   eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
0585   return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
0586 }
0587 
0588 /** \returns the matrix represented by the decomposition,
0589  * i.e., it returns the product: \f$ P^{-1} L U Q^{-1} \f$.
0590  * This function is provided for debug purposes. */
0591 template<typename MatrixType>
0592 MatrixType FullPivLU<MatrixType>::reconstructedMatrix() const
0593 {
0594   eigen_assert(m_isInitialized && "LU is not initialized.");
0595   const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
0596   // LU
0597   MatrixType res(m_lu.rows(),m_lu.cols());
0598   // FIXME the .toDenseMatrix() should not be needed...
0599   res = m_lu.leftCols(smalldim)
0600             .template triangularView<UnitLower>().toDenseMatrix()
0601       * m_lu.topRows(smalldim)
0602             .template triangularView<Upper>().toDenseMatrix();
0603 
0604   // P^{-1}(LU)
0605   res = m_p.inverse() * res;
0606 
0607   // (P^{-1}LU)Q^{-1}
0608   res = res * m_q.inverse();
0609 
0610   return res;
0611 }
0612 
0613 /********* Implementation of kernel() **************************************************/
0614 
0615 namespace internal {
0616 template<typename _MatrixType>
0617 struct kernel_retval<FullPivLU<_MatrixType> >
0618   : kernel_retval_base<FullPivLU<_MatrixType> >
0619 {
0620   EIGEN_MAKE_KERNEL_HELPERS(FullPivLU<_MatrixType>)
0621 
0622   enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
0623             MatrixType::MaxColsAtCompileTime,
0624             MatrixType::MaxRowsAtCompileTime)
0625   };
0626 
0627   template<typename Dest> void evalTo(Dest& dst) const
0628   {
0629     using std::abs;
0630     const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
0631     if(dimker == 0)
0632     {
0633       // The Kernel is just {0}, so it doesn't have a basis properly speaking, but let's
0634       // avoid crashing/asserting as that depends on floating point calculations. Let's
0635       // just return a single column vector filled with zeros.
0636       dst.setZero();
0637       return;
0638     }
0639 
0640     /* Let us use the following lemma:
0641       *
0642       * Lemma: If the matrix A has the LU decomposition PAQ = LU,
0643       * then Ker A = Q(Ker U).
0644       *
0645       * Proof: trivial: just keep in mind that P, Q, L are invertible.
0646       */
0647 
0648     /* Thus, all we need to do is to compute Ker U, and then apply Q.
0649       *
0650       * U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
0651       * Thus, the diagonal of U ends with exactly
0652       * dimKer zero's. Let us use that to construct dimKer linearly
0653       * independent vectors in Ker U.
0654       */
0655 
0656     Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
0657     RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
0658     Index p = 0;
0659     for(Index i = 0; i < dec().nonzeroPivots(); ++i)
0660       if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
0661         pivots.coeffRef(p++) = i;
0662     eigen_internal_assert(p == rank());
0663 
0664     // we construct a temporaty trapezoid matrix m, by taking the U matrix and
0665     // permuting the rows and cols to bring the nonnegligible pivots to the top of
0666     // the main diagonal. We need that to be able to apply our triangular solvers.
0667     // FIXME when we get triangularView-for-rectangular-matrices, this can be simplified
0668     Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
0669            MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
0670       m(dec().matrixLU().block(0, 0, rank(), cols));
0671     for(Index i = 0; i < rank(); ++i)
0672     {
0673       if(i) m.row(i).head(i).setZero();
0674       m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
0675     }
0676     m.block(0, 0, rank(), rank());
0677     m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
0678     for(Index i = 0; i < rank(); ++i)
0679       m.col(i).swap(m.col(pivots.coeff(i)));
0680 
0681     // ok, we have our trapezoid matrix, we can apply the triangular solver.
0682     // notice that the math behind this suggests that we should apply this to the
0683     // negative of the RHS, but for performance we just put the negative sign elsewhere, see below.
0684     m.topLeftCorner(rank(), rank())
0685      .template triangularView<Upper>().solveInPlace(
0686         m.topRightCorner(rank(), dimker)
0687       );
0688 
0689     // now we must undo the column permutation that we had applied!
0690     for(Index i = rank()-1; i >= 0; --i)
0691       m.col(i).swap(m.col(pivots.coeff(i)));
0692 
0693     // see the negative sign in the next line, that's what we were talking about above.
0694     for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
0695     for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
0696     for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
0697   }
0698 };
0699 
0700 /***** Implementation of image() *****************************************************/
0701 
0702 template<typename _MatrixType>
0703 struct image_retval<FullPivLU<_MatrixType> >
0704   : image_retval_base<FullPivLU<_MatrixType> >
0705 {
0706   EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
0707 
0708   enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
0709             MatrixType::MaxColsAtCompileTime,
0710             MatrixType::MaxRowsAtCompileTime)
0711   };
0712 
0713   template<typename Dest> void evalTo(Dest& dst) const
0714   {
0715     using std::abs;
0716     if(rank() == 0)
0717     {
0718       // The Image is just {0}, so it doesn't have a basis properly speaking, but let's
0719       // avoid crashing/asserting as that depends on floating point calculations. Let's
0720       // just return a single column vector filled with zeros.
0721       dst.setZero();
0722       return;
0723     }
0724 
0725     Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
0726     RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
0727     Index p = 0;
0728     for(Index i = 0; i < dec().nonzeroPivots(); ++i)
0729       if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
0730         pivots.coeffRef(p++) = i;
0731     eigen_internal_assert(p == rank());
0732 
0733     for(Index i = 0; i < rank(); ++i)
0734       dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
0735   }
0736 };
0737 
0738 /***** Implementation of solve() *****************************************************/
0739 
0740 } // end namespace internal
0741 
0742 #ifndef EIGEN_PARSED_BY_DOXYGEN
0743 template<typename _MatrixType>
0744 template<typename RhsType, typename DstType>
0745 void FullPivLU<_MatrixType>::_solve_impl(const RhsType &rhs, DstType &dst) const
0746 {
0747   /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
0748   * So we proceed as follows:
0749   * Step 1: compute c = P * rhs.
0750   * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
0751   * Step 3: replace c by the solution x to Ux = c. May or may not exist.
0752   * Step 4: result = Q * c;
0753   */
0754 
0755   const Index rows = this->rows(),
0756               cols = this->cols(),
0757               nonzero_pivots = this->rank();
0758   const Index smalldim = (std::min)(rows, cols);
0759 
0760   if(nonzero_pivots == 0)
0761   {
0762     dst.setZero();
0763     return;
0764   }
0765 
0766   typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
0767 
0768   // Step 1
0769   c = permutationP() * rhs;
0770 
0771   // Step 2
0772   m_lu.topLeftCorner(smalldim,smalldim)
0773       .template triangularView<UnitLower>()
0774       .solveInPlace(c.topRows(smalldim));
0775   if(rows>cols)
0776     c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
0777 
0778   // Step 3
0779   m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
0780       .template triangularView<Upper>()
0781       .solveInPlace(c.topRows(nonzero_pivots));
0782 
0783   // Step 4
0784   for(Index i = 0; i < nonzero_pivots; ++i)
0785     dst.row(permutationQ().indices().coeff(i)) = c.row(i);
0786   for(Index i = nonzero_pivots; i < m_lu.cols(); ++i)
0787     dst.row(permutationQ().indices().coeff(i)).setZero();
0788 }
0789 
0790 template<typename _MatrixType>
0791 template<bool Conjugate, typename RhsType, typename DstType>
0792 void FullPivLU<_MatrixType>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
0793 {
0794   /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1},
0795    * and since permutations are real and unitary, we can write this
0796    * as   A^T = Q U^T L^T P,
0797    * So we proceed as follows:
0798    * Step 1: compute c = Q^T rhs.
0799    * Step 2: replace c by the solution x to U^T x = c. May or may not exist.
0800    * Step 3: replace c by the solution x to L^T x = c.
0801    * Step 4: result = P^T c.
0802    * If Conjugate is true, replace "^T" by "^*" above.
0803    */
0804 
0805   const Index rows = this->rows(), cols = this->cols(),
0806     nonzero_pivots = this->rank();
0807   const Index smalldim = (std::min)(rows, cols);
0808 
0809   if(nonzero_pivots == 0)
0810   {
0811     dst.setZero();
0812     return;
0813   }
0814 
0815   typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
0816 
0817   // Step 1
0818   c = permutationQ().inverse() * rhs;
0819 
0820   // Step 2
0821   m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
0822       .template triangularView<Upper>()
0823       .transpose()
0824       .template conjugateIf<Conjugate>()
0825       .solveInPlace(c.topRows(nonzero_pivots));
0826 
0827   // Step 3
0828   m_lu.topLeftCorner(smalldim, smalldim)
0829       .template triangularView<UnitLower>()
0830       .transpose()
0831       .template conjugateIf<Conjugate>()
0832       .solveInPlace(c.topRows(smalldim));
0833 
0834   // Step 4
0835   PermutationPType invp = permutationP().inverse().eval();
0836   for(Index i = 0; i < smalldim; ++i)
0837     dst.row(invp.indices().coeff(i)) = c.row(i);
0838   for(Index i = smalldim; i < rows; ++i)
0839     dst.row(invp.indices().coeff(i)).setZero();
0840 }
0841 
0842 #endif
0843 
0844 namespace internal {
0845 
0846 
0847 /***** Implementation of inverse() *****************************************************/
0848 template<typename DstXprType, typename MatrixType>
0849 struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename FullPivLU<MatrixType>::Scalar>, Dense2Dense>
0850 {
0851   typedef FullPivLU<MatrixType> LuType;
0852   typedef Inverse<LuType> SrcXprType;
0853   static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename MatrixType::Scalar> &)
0854   {
0855     dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
0856   }
0857 };
0858 } // end namespace internal
0859 
0860 /******* MatrixBase methods *****************************************************************/
0861 
0862 /** \lu_module
0863   *
0864   * \return the full-pivoting LU decomposition of \c *this.
0865   *
0866   * \sa class FullPivLU
0867   */
0868 template<typename Derived>
0869 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
0870 MatrixBase<Derived>::fullPivLu() const
0871 {
0872   return FullPivLU<PlainObject>(eval());
0873 }
0874 
0875 } // end namespace Eigen
0876 
0877 #endif // EIGEN_LU_H