1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library 2*bf2c3715SXin Li // for linear algebra. 3*bf2c3715SXin Li // 4*bf2c3715SXin Li // Copyright (C) 2008-2011 Gael Guennebaud <[email protected]> 5*bf2c3715SXin Li // Copyright (C) 2009 Keir Mierle <[email protected]> 6*bf2c3715SXin Li // Copyright (C) 2009 Benoit Jacob <[email protected]> 7*bf2c3715SXin Li // Copyright (C) 2011 Timothy E. Holy <[email protected] > 8*bf2c3715SXin Li // 9*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla 10*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed 11*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 12*bf2c3715SXin Li 13*bf2c3715SXin Li #ifndef EIGEN_LDLT_H 14*bf2c3715SXin Li #define EIGEN_LDLT_H 15*bf2c3715SXin Li 16*bf2c3715SXin Li namespace Eigen { 17*bf2c3715SXin Li 18*bf2c3715SXin Li namespace internal { 19*bf2c3715SXin Li template<typename _MatrixType, int _UpLo> struct traits<LDLT<_MatrixType, _UpLo> > 20*bf2c3715SXin Li : traits<_MatrixType> 21*bf2c3715SXin Li { 22*bf2c3715SXin Li typedef MatrixXpr XprKind; 23*bf2c3715SXin Li typedef SolverStorage StorageKind; 24*bf2c3715SXin Li typedef int StorageIndex; 25*bf2c3715SXin Li enum { Flags = 0 }; 26*bf2c3715SXin Li }; 27*bf2c3715SXin Li 28*bf2c3715SXin Li template<typename MatrixType, int UpLo> struct LDLT_Traits; 29*bf2c3715SXin Li 30*bf2c3715SXin Li // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef 31*bf2c3715SXin Li enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite }; 32*bf2c3715SXin Li } 33*bf2c3715SXin Li 34*bf2c3715SXin Li /** \ingroup Cholesky_Module 35*bf2c3715SXin Li * 36*bf2c3715SXin Li * \class LDLT 37*bf2c3715SXin Li * 38*bf2c3715SXin Li * \brief Robust Cholesky decomposition of a matrix with pivoting 39*bf2c3715SXin Li * 40*bf2c3715SXin Li * \tparam _MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition 41*bf2c3715SXin Li * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper. 42*bf2c3715SXin Li * The other triangular part won't be read. 43*bf2c3715SXin Li * 44*bf2c3715SXin Li * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite 45*bf2c3715SXin Li * matrix \f$ A \f$ such that \f$ A = P^TLDL^*P \f$, where P is a permutation matrix, L 46*bf2c3715SXin Li * is lower triangular with a unit diagonal and D is a diagonal matrix. 47*bf2c3715SXin Li * 48*bf2c3715SXin Li * The decomposition uses pivoting to ensure stability, so that D will have 49*bf2c3715SXin Li * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root 50*bf2c3715SXin Li * on D also stabilizes the computation. 51*bf2c3715SXin Li * 52*bf2c3715SXin Li * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky 53*bf2c3715SXin Li * decomposition to determine whether a system of equations has a solution. 54*bf2c3715SXin Li * 55*bf2c3715SXin Li * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism. 56*bf2c3715SXin Li * 57*bf2c3715SXin Li * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT 58*bf2c3715SXin Li */ 59*bf2c3715SXin Li template<typename _MatrixType, int _UpLo> class LDLT 60*bf2c3715SXin Li : public SolverBase<LDLT<_MatrixType, _UpLo> > 61*bf2c3715SXin Li { 62*bf2c3715SXin Li public: 63*bf2c3715SXin Li typedef _MatrixType MatrixType; 64*bf2c3715SXin Li typedef SolverBase<LDLT> Base; 65*bf2c3715SXin Li friend class SolverBase<LDLT>; 66*bf2c3715SXin Li 67*bf2c3715SXin Li EIGEN_GENERIC_PUBLIC_INTERFACE(LDLT) 68*bf2c3715SXin Li enum { 69*bf2c3715SXin Li MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, 70*bf2c3715SXin Li MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, 71*bf2c3715SXin Li UpLo = _UpLo 72*bf2c3715SXin Li }; 73*bf2c3715SXin Li typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> TmpMatrixType; 74*bf2c3715SXin Li 75*bf2c3715SXin Li typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType; 76*bf2c3715SXin Li typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType; 77*bf2c3715SXin Li 78*bf2c3715SXin Li typedef internal::LDLT_Traits<MatrixType,UpLo> Traits; 79*bf2c3715SXin Li 80*bf2c3715SXin Li /** \brief Default Constructor. 81*bf2c3715SXin Li * 82*bf2c3715SXin Li * The default constructor is useful in cases in which the user intends to 83*bf2c3715SXin Li * perform decompositions via LDLT::compute(const MatrixType&). 84*bf2c3715SXin Li */ 85*bf2c3715SXin Li LDLT() 86*bf2c3715SXin Li : m_matrix(), 87*bf2c3715SXin Li m_transpositions(), 88*bf2c3715SXin Li m_sign(internal::ZeroSign), 89*bf2c3715SXin Li m_isInitialized(false) 90*bf2c3715SXin Li {} 91*bf2c3715SXin Li 92*bf2c3715SXin Li /** \brief Default Constructor with memory preallocation 93*bf2c3715SXin Li * 94*bf2c3715SXin Li * Like the default constructor but with preallocation of the internal data 95*bf2c3715SXin Li * according to the specified problem \a size. 96*bf2c3715SXin Li * \sa LDLT() 97*bf2c3715SXin Li */ 98*bf2c3715SXin Li explicit LDLT(Index size) 99*bf2c3715SXin Li : m_matrix(size, size), 100*bf2c3715SXin Li m_transpositions(size), 101*bf2c3715SXin Li m_temporary(size), 102*bf2c3715SXin Li m_sign(internal::ZeroSign), 103*bf2c3715SXin Li m_isInitialized(false) 104*bf2c3715SXin Li {} 105*bf2c3715SXin Li 106*bf2c3715SXin Li /** \brief Constructor with decomposition 107*bf2c3715SXin Li * 108*bf2c3715SXin Li * This calculates the decomposition for the input \a matrix. 109*bf2c3715SXin Li * 110*bf2c3715SXin Li * \sa LDLT(Index size) 111*bf2c3715SXin Li */ 112*bf2c3715SXin Li template<typename InputType> 113*bf2c3715SXin Li explicit LDLT(const EigenBase<InputType>& matrix) 114*bf2c3715SXin Li : m_matrix(matrix.rows(), matrix.cols()), 115*bf2c3715SXin Li m_transpositions(matrix.rows()), 116*bf2c3715SXin Li m_temporary(matrix.rows()), 117*bf2c3715SXin Li m_sign(internal::ZeroSign), 118*bf2c3715SXin Li m_isInitialized(false) 119*bf2c3715SXin Li { 120*bf2c3715SXin Li compute(matrix.derived()); 121*bf2c3715SXin Li } 122*bf2c3715SXin Li 123*bf2c3715SXin Li /** \brief Constructs a LDLT factorization from a given matrix 124*bf2c3715SXin Li * 125*bf2c3715SXin Li * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref. 126*bf2c3715SXin Li * 127*bf2c3715SXin Li * \sa LDLT(const EigenBase&) 128*bf2c3715SXin Li */ 129*bf2c3715SXin Li template<typename InputType> 130*bf2c3715SXin Li explicit LDLT(EigenBase<InputType>& matrix) 131*bf2c3715SXin Li : m_matrix(matrix.derived()), 132*bf2c3715SXin Li m_transpositions(matrix.rows()), 133*bf2c3715SXin Li m_temporary(matrix.rows()), 134*bf2c3715SXin Li m_sign(internal::ZeroSign), 135*bf2c3715SXin Li m_isInitialized(false) 136*bf2c3715SXin Li { 137*bf2c3715SXin Li compute(matrix.derived()); 138*bf2c3715SXin Li } 139*bf2c3715SXin Li 140*bf2c3715SXin Li /** Clear any existing decomposition 141*bf2c3715SXin Li * \sa rankUpdate(w,sigma) 142*bf2c3715SXin Li */ 143*bf2c3715SXin Li void setZero() 144*bf2c3715SXin Li { 145*bf2c3715SXin Li m_isInitialized = false; 146*bf2c3715SXin Li } 147*bf2c3715SXin Li 148*bf2c3715SXin Li /** \returns a view of the upper triangular matrix U */ 149*bf2c3715SXin Li inline typename Traits::MatrixU matrixU() const 150*bf2c3715SXin Li { 151*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 152*bf2c3715SXin Li return Traits::getU(m_matrix); 153*bf2c3715SXin Li } 154*bf2c3715SXin Li 155*bf2c3715SXin Li /** \returns a view of the lower triangular matrix L */ 156*bf2c3715SXin Li inline typename Traits::MatrixL matrixL() const 157*bf2c3715SXin Li { 158*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 159*bf2c3715SXin Li return Traits::getL(m_matrix); 160*bf2c3715SXin Li } 161*bf2c3715SXin Li 162*bf2c3715SXin Li /** \returns the permutation matrix P as a transposition sequence. 163*bf2c3715SXin Li */ 164*bf2c3715SXin Li inline const TranspositionType& transpositionsP() const 165*bf2c3715SXin Li { 166*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 167*bf2c3715SXin Li return m_transpositions; 168*bf2c3715SXin Li } 169*bf2c3715SXin Li 170*bf2c3715SXin Li /** \returns the coefficients of the diagonal matrix D */ 171*bf2c3715SXin Li inline Diagonal<const MatrixType> vectorD() const 172*bf2c3715SXin Li { 173*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 174*bf2c3715SXin Li return m_matrix.diagonal(); 175*bf2c3715SXin Li } 176*bf2c3715SXin Li 177*bf2c3715SXin Li /** \returns true if the matrix is positive (semidefinite) */ 178*bf2c3715SXin Li inline bool isPositive() const 179*bf2c3715SXin Li { 180*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 181*bf2c3715SXin Li return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign; 182*bf2c3715SXin Li } 183*bf2c3715SXin Li 184*bf2c3715SXin Li /** \returns true if the matrix is negative (semidefinite) */ 185*bf2c3715SXin Li inline bool isNegative(void) const 186*bf2c3715SXin Li { 187*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 188*bf2c3715SXin Li return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign; 189*bf2c3715SXin Li } 190*bf2c3715SXin Li 191*bf2c3715SXin Li #ifdef EIGEN_PARSED_BY_DOXYGEN 192*bf2c3715SXin Li /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A. 193*bf2c3715SXin Li * 194*bf2c3715SXin Li * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> . 195*bf2c3715SXin Li * 196*bf2c3715SXin Li * \note_about_checking_solutions 197*bf2c3715SXin Li * 198*bf2c3715SXin Li * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$ 199*bf2c3715SXin Li * by solving the systems \f$ P^T y_1 = b \f$, \f$ L y_2 = y_1 \f$, \f$ D y_3 = y_2 \f$, 200*bf2c3715SXin Li * \f$ L^* y_4 = y_3 \f$ and \f$ P x = y_4 \f$ in succession. If the matrix \f$ A \f$ is singular, then 201*bf2c3715SXin Li * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the 202*bf2c3715SXin Li * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function 203*bf2c3715SXin Li * computes the least-square solution of \f$ A x = b \f$ if \f$ A \f$ is singular. 204*bf2c3715SXin Li * 205*bf2c3715SXin Li * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt() 206*bf2c3715SXin Li */ 207*bf2c3715SXin Li template<typename Rhs> 208*bf2c3715SXin Li inline const Solve<LDLT, Rhs> 209*bf2c3715SXin Li solve(const MatrixBase<Rhs>& b) const; 210*bf2c3715SXin Li #endif 211*bf2c3715SXin Li 212*bf2c3715SXin Li template<typename Derived> 213*bf2c3715SXin Li bool solveInPlace(MatrixBase<Derived> &bAndX) const; 214*bf2c3715SXin Li 215*bf2c3715SXin Li template<typename InputType> 216*bf2c3715SXin Li LDLT& compute(const EigenBase<InputType>& matrix); 217*bf2c3715SXin Li 218*bf2c3715SXin Li /** \returns an estimate of the reciprocal condition number of the matrix of 219*bf2c3715SXin Li * which \c *this is the LDLT decomposition. 220*bf2c3715SXin Li */ 221*bf2c3715SXin Li RealScalar rcond() const 222*bf2c3715SXin Li { 223*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 224*bf2c3715SXin Li return internal::rcond_estimate_helper(m_l1_norm, *this); 225*bf2c3715SXin Li } 226*bf2c3715SXin Li 227*bf2c3715SXin Li template <typename Derived> 228*bf2c3715SXin Li LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1); 229*bf2c3715SXin Li 230*bf2c3715SXin Li /** \returns the internal LDLT decomposition matrix 231*bf2c3715SXin Li * 232*bf2c3715SXin Li * TODO: document the storage layout 233*bf2c3715SXin Li */ 234*bf2c3715SXin Li inline const MatrixType& matrixLDLT() const 235*bf2c3715SXin Li { 236*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 237*bf2c3715SXin Li return m_matrix; 238*bf2c3715SXin Li } 239*bf2c3715SXin Li 240*bf2c3715SXin Li MatrixType reconstructedMatrix() const; 241*bf2c3715SXin Li 242*bf2c3715SXin Li /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint. 243*bf2c3715SXin Li * 244*bf2c3715SXin Li * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as: 245*bf2c3715SXin Li * \code x = decomposition.adjoint().solve(b) \endcode 246*bf2c3715SXin Li */ 247*bf2c3715SXin Li const LDLT& adjoint() const { return *this; }; 248*bf2c3715SXin Li 249*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); } 250*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); } 251*bf2c3715SXin Li 252*bf2c3715SXin Li /** \brief Reports whether previous computation was successful. 253*bf2c3715SXin Li * 254*bf2c3715SXin Li * \returns \c Success if computation was successful, 255*bf2c3715SXin Li * \c NumericalIssue if the factorization failed because of a zero pivot. 256*bf2c3715SXin Li */ 257*bf2c3715SXin Li ComputationInfo info() const 258*bf2c3715SXin Li { 259*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 260*bf2c3715SXin Li return m_info; 261*bf2c3715SXin Li } 262*bf2c3715SXin Li 263*bf2c3715SXin Li #ifndef EIGEN_PARSED_BY_DOXYGEN 264*bf2c3715SXin Li template<typename RhsType, typename DstType> 265*bf2c3715SXin Li void _solve_impl(const RhsType &rhs, DstType &dst) const; 266*bf2c3715SXin Li 267*bf2c3715SXin Li template<bool Conjugate, typename RhsType, typename DstType> 268*bf2c3715SXin Li void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const; 269*bf2c3715SXin Li #endif 270*bf2c3715SXin Li 271*bf2c3715SXin Li protected: 272*bf2c3715SXin Li 273*bf2c3715SXin Li static void check_template_parameters() 274*bf2c3715SXin Li { 275*bf2c3715SXin Li EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar); 276*bf2c3715SXin Li } 277*bf2c3715SXin Li 278*bf2c3715SXin Li /** \internal 279*bf2c3715SXin Li * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U. 280*bf2c3715SXin Li * The strict upper part is used during the decomposition, the strict lower 281*bf2c3715SXin Li * part correspond to the coefficients of L (its diagonal is equal to 1 and 282*bf2c3715SXin Li * is not stored), and the diagonal entries correspond to D. 283*bf2c3715SXin Li */ 284*bf2c3715SXin Li MatrixType m_matrix; 285*bf2c3715SXin Li RealScalar m_l1_norm; 286*bf2c3715SXin Li TranspositionType m_transpositions; 287*bf2c3715SXin Li TmpMatrixType m_temporary; 288*bf2c3715SXin Li internal::SignMatrix m_sign; 289*bf2c3715SXin Li bool m_isInitialized; 290*bf2c3715SXin Li ComputationInfo m_info; 291*bf2c3715SXin Li }; 292*bf2c3715SXin Li 293*bf2c3715SXin Li namespace internal { 294*bf2c3715SXin Li 295*bf2c3715SXin Li template<int UpLo> struct ldlt_inplace; 296*bf2c3715SXin Li 297*bf2c3715SXin Li template<> struct ldlt_inplace<Lower> 298*bf2c3715SXin Li { 299*bf2c3715SXin Li template<typename MatrixType, typename TranspositionType, typename Workspace> 300*bf2c3715SXin Li static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign) 301*bf2c3715SXin Li { 302*bf2c3715SXin Li using std::abs; 303*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar; 304*bf2c3715SXin Li typedef typename MatrixType::RealScalar RealScalar; 305*bf2c3715SXin Li typedef typename TranspositionType::StorageIndex IndexType; 306*bf2c3715SXin Li eigen_assert(mat.rows()==mat.cols()); 307*bf2c3715SXin Li const Index size = mat.rows(); 308*bf2c3715SXin Li bool found_zero_pivot = false; 309*bf2c3715SXin Li bool ret = true; 310*bf2c3715SXin Li 311*bf2c3715SXin Li if (size <= 1) 312*bf2c3715SXin Li { 313*bf2c3715SXin Li transpositions.setIdentity(); 314*bf2c3715SXin Li if(size==0) sign = ZeroSign; 315*bf2c3715SXin Li else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef; 316*bf2c3715SXin Li else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef; 317*bf2c3715SXin Li else sign = ZeroSign; 318*bf2c3715SXin Li return true; 319*bf2c3715SXin Li } 320*bf2c3715SXin Li 321*bf2c3715SXin Li for (Index k = 0; k < size; ++k) 322*bf2c3715SXin Li { 323*bf2c3715SXin Li // Find largest diagonal element 324*bf2c3715SXin Li Index index_of_biggest_in_corner; 325*bf2c3715SXin Li mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner); 326*bf2c3715SXin Li index_of_biggest_in_corner += k; 327*bf2c3715SXin Li 328*bf2c3715SXin Li transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner); 329*bf2c3715SXin Li if(k != index_of_biggest_in_corner) 330*bf2c3715SXin Li { 331*bf2c3715SXin Li // apply the transposition while taking care to consider only 332*bf2c3715SXin Li // the lower triangular part 333*bf2c3715SXin Li Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element 334*bf2c3715SXin Li mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k)); 335*bf2c3715SXin Li mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s)); 336*bf2c3715SXin Li std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner)); 337*bf2c3715SXin Li for(Index i=k+1;i<index_of_biggest_in_corner;++i) 338*bf2c3715SXin Li { 339*bf2c3715SXin Li Scalar tmp = mat.coeffRef(i,k); 340*bf2c3715SXin Li mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i)); 341*bf2c3715SXin Li mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp); 342*bf2c3715SXin Li } 343*bf2c3715SXin Li if(NumTraits<Scalar>::IsComplex) 344*bf2c3715SXin Li mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k)); 345*bf2c3715SXin Li } 346*bf2c3715SXin Li 347*bf2c3715SXin Li // partition the matrix: 348*bf2c3715SXin Li // A00 | - | - 349*bf2c3715SXin Li // lu = A10 | A11 | - 350*bf2c3715SXin Li // A20 | A21 | A22 351*bf2c3715SXin Li Index rs = size - k - 1; 352*bf2c3715SXin Li Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1); 353*bf2c3715SXin Li Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k); 354*bf2c3715SXin Li Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k); 355*bf2c3715SXin Li 356*bf2c3715SXin Li if(k>0) 357*bf2c3715SXin Li { 358*bf2c3715SXin Li temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint(); 359*bf2c3715SXin Li mat.coeffRef(k,k) -= (A10 * temp.head(k)).value(); 360*bf2c3715SXin Li if(rs>0) 361*bf2c3715SXin Li A21.noalias() -= A20 * temp.head(k); 362*bf2c3715SXin Li } 363*bf2c3715SXin Li 364*bf2c3715SXin Li // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot 365*bf2c3715SXin Li // was smaller than the cutoff value. However, since LDLT is not rank-revealing 366*bf2c3715SXin Li // we should only make sure that we do not introduce INF or NaN values. 367*bf2c3715SXin Li // Remark that LAPACK also uses 0 as the cutoff value. 368*bf2c3715SXin Li RealScalar realAkk = numext::real(mat.coeffRef(k,k)); 369*bf2c3715SXin Li bool pivot_is_valid = (abs(realAkk) > RealScalar(0)); 370*bf2c3715SXin Li 371*bf2c3715SXin Li if(k==0 && !pivot_is_valid) 372*bf2c3715SXin Li { 373*bf2c3715SXin Li // The entire diagonal is zero, there is nothing more to do 374*bf2c3715SXin Li // except filling the transpositions, and checking whether the matrix is zero. 375*bf2c3715SXin Li sign = ZeroSign; 376*bf2c3715SXin Li for(Index j = 0; j<size; ++j) 377*bf2c3715SXin Li { 378*bf2c3715SXin Li transpositions.coeffRef(j) = IndexType(j); 379*bf2c3715SXin Li ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all(); 380*bf2c3715SXin Li } 381*bf2c3715SXin Li return ret; 382*bf2c3715SXin Li } 383*bf2c3715SXin Li 384*bf2c3715SXin Li if((rs>0) && pivot_is_valid) 385*bf2c3715SXin Li A21 /= realAkk; 386*bf2c3715SXin Li else if(rs>0) 387*bf2c3715SXin Li ret = ret && (A21.array()==Scalar(0)).all(); 388*bf2c3715SXin Li 389*bf2c3715SXin Li if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed 390*bf2c3715SXin Li else if(!pivot_is_valid) found_zero_pivot = true; 391*bf2c3715SXin Li 392*bf2c3715SXin Li if (sign == PositiveSemiDef) { 393*bf2c3715SXin Li if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite; 394*bf2c3715SXin Li } else if (sign == NegativeSemiDef) { 395*bf2c3715SXin Li if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite; 396*bf2c3715SXin Li } else if (sign == ZeroSign) { 397*bf2c3715SXin Li if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef; 398*bf2c3715SXin Li else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef; 399*bf2c3715SXin Li } 400*bf2c3715SXin Li } 401*bf2c3715SXin Li 402*bf2c3715SXin Li return ret; 403*bf2c3715SXin Li } 404*bf2c3715SXin Li 405*bf2c3715SXin Li // Reference for the algorithm: Davis and Hager, "Multiple Rank 406*bf2c3715SXin Li // Modifications of a Sparse Cholesky Factorization" (Algorithm 1) 407*bf2c3715SXin Li // Trivial rearrangements of their computations (Timothy E. Holy) 408*bf2c3715SXin Li // allow their algorithm to work for rank-1 updates even if the 409*bf2c3715SXin Li // original matrix is not of full rank. 410*bf2c3715SXin Li // Here only rank-1 updates are implemented, to reduce the 411*bf2c3715SXin Li // requirement for intermediate storage and improve accuracy 412*bf2c3715SXin Li template<typename MatrixType, typename WDerived> 413*bf2c3715SXin Li static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1) 414*bf2c3715SXin Li { 415*bf2c3715SXin Li using numext::isfinite; 416*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar; 417*bf2c3715SXin Li typedef typename MatrixType::RealScalar RealScalar; 418*bf2c3715SXin Li 419*bf2c3715SXin Li const Index size = mat.rows(); 420*bf2c3715SXin Li eigen_assert(mat.cols() == size && w.size()==size); 421*bf2c3715SXin Li 422*bf2c3715SXin Li RealScalar alpha = 1; 423*bf2c3715SXin Li 424*bf2c3715SXin Li // Apply the update 425*bf2c3715SXin Li for (Index j = 0; j < size; j++) 426*bf2c3715SXin Li { 427*bf2c3715SXin Li // Check for termination due to an original decomposition of low-rank 428*bf2c3715SXin Li if (!(isfinite)(alpha)) 429*bf2c3715SXin Li break; 430*bf2c3715SXin Li 431*bf2c3715SXin Li // Update the diagonal terms 432*bf2c3715SXin Li RealScalar dj = numext::real(mat.coeff(j,j)); 433*bf2c3715SXin Li Scalar wj = w.coeff(j); 434*bf2c3715SXin Li RealScalar swj2 = sigma*numext::abs2(wj); 435*bf2c3715SXin Li RealScalar gamma = dj*alpha + swj2; 436*bf2c3715SXin Li 437*bf2c3715SXin Li mat.coeffRef(j,j) += swj2/alpha; 438*bf2c3715SXin Li alpha += swj2/dj; 439*bf2c3715SXin Li 440*bf2c3715SXin Li 441*bf2c3715SXin Li // Update the terms of L 442*bf2c3715SXin Li Index rs = size-j-1; 443*bf2c3715SXin Li w.tail(rs) -= wj * mat.col(j).tail(rs); 444*bf2c3715SXin Li if(gamma != 0) 445*bf2c3715SXin Li mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs); 446*bf2c3715SXin Li } 447*bf2c3715SXin Li return true; 448*bf2c3715SXin Li } 449*bf2c3715SXin Li 450*bf2c3715SXin Li template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType> 451*bf2c3715SXin Li static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1) 452*bf2c3715SXin Li { 453*bf2c3715SXin Li // Apply the permutation to the input w 454*bf2c3715SXin Li tmp = transpositions * w; 455*bf2c3715SXin Li 456*bf2c3715SXin Li return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma); 457*bf2c3715SXin Li } 458*bf2c3715SXin Li }; 459*bf2c3715SXin Li 460*bf2c3715SXin Li template<> struct ldlt_inplace<Upper> 461*bf2c3715SXin Li { 462*bf2c3715SXin Li template<typename MatrixType, typename TranspositionType, typename Workspace> 463*bf2c3715SXin Li static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign) 464*bf2c3715SXin Li { 465*bf2c3715SXin Li Transpose<MatrixType> matt(mat); 466*bf2c3715SXin Li return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign); 467*bf2c3715SXin Li } 468*bf2c3715SXin Li 469*bf2c3715SXin Li template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType> 470*bf2c3715SXin Li static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1) 471*bf2c3715SXin Li { 472*bf2c3715SXin Li Transpose<MatrixType> matt(mat); 473*bf2c3715SXin Li return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma); 474*bf2c3715SXin Li } 475*bf2c3715SXin Li }; 476*bf2c3715SXin Li 477*bf2c3715SXin Li template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower> 478*bf2c3715SXin Li { 479*bf2c3715SXin Li typedef const TriangularView<const MatrixType, UnitLower> MatrixL; 480*bf2c3715SXin Li typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU; 481*bf2c3715SXin Li static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); } 482*bf2c3715SXin Li static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); } 483*bf2c3715SXin Li }; 484*bf2c3715SXin Li 485*bf2c3715SXin Li template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper> 486*bf2c3715SXin Li { 487*bf2c3715SXin Li typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL; 488*bf2c3715SXin Li typedef const TriangularView<const MatrixType, UnitUpper> MatrixU; 489*bf2c3715SXin Li static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); } 490*bf2c3715SXin Li static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); } 491*bf2c3715SXin Li }; 492*bf2c3715SXin Li 493*bf2c3715SXin Li } // end namespace internal 494*bf2c3715SXin Li 495*bf2c3715SXin Li /** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix 496*bf2c3715SXin Li */ 497*bf2c3715SXin Li template<typename MatrixType, int _UpLo> 498*bf2c3715SXin Li template<typename InputType> 499*bf2c3715SXin Li LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a) 500*bf2c3715SXin Li { 501*bf2c3715SXin Li check_template_parameters(); 502*bf2c3715SXin Li 503*bf2c3715SXin Li eigen_assert(a.rows()==a.cols()); 504*bf2c3715SXin Li const Index size = a.rows(); 505*bf2c3715SXin Li 506*bf2c3715SXin Li m_matrix = a.derived(); 507*bf2c3715SXin Li 508*bf2c3715SXin Li // Compute matrix L1 norm = max abs column sum. 509*bf2c3715SXin Li m_l1_norm = RealScalar(0); 510*bf2c3715SXin Li // TODO move this code to SelfAdjointView 511*bf2c3715SXin Li for (Index col = 0; col < size; ++col) { 512*bf2c3715SXin Li RealScalar abs_col_sum; 513*bf2c3715SXin Li if (_UpLo == Lower) 514*bf2c3715SXin Li abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>(); 515*bf2c3715SXin Li else 516*bf2c3715SXin Li abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>(); 517*bf2c3715SXin Li if (abs_col_sum > m_l1_norm) 518*bf2c3715SXin Li m_l1_norm = abs_col_sum; 519*bf2c3715SXin Li } 520*bf2c3715SXin Li 521*bf2c3715SXin Li m_transpositions.resize(size); 522*bf2c3715SXin Li m_isInitialized = false; 523*bf2c3715SXin Li m_temporary.resize(size); 524*bf2c3715SXin Li m_sign = internal::ZeroSign; 525*bf2c3715SXin Li 526*bf2c3715SXin Li m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue; 527*bf2c3715SXin Li 528*bf2c3715SXin Li m_isInitialized = true; 529*bf2c3715SXin Li return *this; 530*bf2c3715SXin Li } 531*bf2c3715SXin Li 532*bf2c3715SXin Li /** Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T. 533*bf2c3715SXin Li * \param w a vector to be incorporated into the decomposition. 534*bf2c3715SXin Li * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1. 535*bf2c3715SXin Li * \sa setZero() 536*bf2c3715SXin Li */ 537*bf2c3715SXin Li template<typename MatrixType, int _UpLo> 538*bf2c3715SXin Li template<typename Derived> 539*bf2c3715SXin Li LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma) 540*bf2c3715SXin Li { 541*bf2c3715SXin Li typedef typename TranspositionType::StorageIndex IndexType; 542*bf2c3715SXin Li const Index size = w.rows(); 543*bf2c3715SXin Li if (m_isInitialized) 544*bf2c3715SXin Li { 545*bf2c3715SXin Li eigen_assert(m_matrix.rows()==size); 546*bf2c3715SXin Li } 547*bf2c3715SXin Li else 548*bf2c3715SXin Li { 549*bf2c3715SXin Li m_matrix.resize(size,size); 550*bf2c3715SXin Li m_matrix.setZero(); 551*bf2c3715SXin Li m_transpositions.resize(size); 552*bf2c3715SXin Li for (Index i = 0; i < size; i++) 553*bf2c3715SXin Li m_transpositions.coeffRef(i) = IndexType(i); 554*bf2c3715SXin Li m_temporary.resize(size); 555*bf2c3715SXin Li m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef; 556*bf2c3715SXin Li m_isInitialized = true; 557*bf2c3715SXin Li } 558*bf2c3715SXin Li 559*bf2c3715SXin Li internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma); 560*bf2c3715SXin Li 561*bf2c3715SXin Li return *this; 562*bf2c3715SXin Li } 563*bf2c3715SXin Li 564*bf2c3715SXin Li #ifndef EIGEN_PARSED_BY_DOXYGEN 565*bf2c3715SXin Li template<typename _MatrixType, int _UpLo> 566*bf2c3715SXin Li template<typename RhsType, typename DstType> 567*bf2c3715SXin Li void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const 568*bf2c3715SXin Li { 569*bf2c3715SXin Li _solve_impl_transposed<true>(rhs, dst); 570*bf2c3715SXin Li } 571*bf2c3715SXin Li 572*bf2c3715SXin Li template<typename _MatrixType,int _UpLo> 573*bf2c3715SXin Li template<bool Conjugate, typename RhsType, typename DstType> 574*bf2c3715SXin Li void LDLT<_MatrixType,_UpLo>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const 575*bf2c3715SXin Li { 576*bf2c3715SXin Li // dst = P b 577*bf2c3715SXin Li dst = m_transpositions * rhs; 578*bf2c3715SXin Li 579*bf2c3715SXin Li // dst = L^-1 (P b) 580*bf2c3715SXin Li // dst = L^-*T (P b) 581*bf2c3715SXin Li matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst); 582*bf2c3715SXin Li 583*bf2c3715SXin Li // dst = D^-* (L^-1 P b) 584*bf2c3715SXin Li // dst = D^-1 (L^-*T P b) 585*bf2c3715SXin Li // more precisely, use pseudo-inverse of D (see bug 241) 586*bf2c3715SXin Li using std::abs; 587*bf2c3715SXin Li const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD()); 588*bf2c3715SXin Li // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min()) 589*bf2c3715SXin Li // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS: 590*bf2c3715SXin Li // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest()); 591*bf2c3715SXin Li // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest 592*bf2c3715SXin Li // diagonal element is not well justified and leads to numerical issues in some cases. 593*bf2c3715SXin Li // Moreover, Lapack's xSYTRS routines use 0 for the tolerance. 594*bf2c3715SXin Li // Using numeric_limits::min() gives us more robustness to denormals. 595*bf2c3715SXin Li RealScalar tolerance = (std::numeric_limits<RealScalar>::min)(); 596*bf2c3715SXin Li for (Index i = 0; i < vecD.size(); ++i) 597*bf2c3715SXin Li { 598*bf2c3715SXin Li if(abs(vecD(i)) > tolerance) 599*bf2c3715SXin Li dst.row(i) /= vecD(i); 600*bf2c3715SXin Li else 601*bf2c3715SXin Li dst.row(i).setZero(); 602*bf2c3715SXin Li } 603*bf2c3715SXin Li 604*bf2c3715SXin Li // dst = L^-* (D^-* L^-1 P b) 605*bf2c3715SXin Li // dst = L^-T (D^-1 L^-*T P b) 606*bf2c3715SXin Li matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst); 607*bf2c3715SXin Li 608*bf2c3715SXin Li // dst = P^T (L^-* D^-* L^-1 P b) = A^-1 b 609*bf2c3715SXin Li // dst = P^-T (L^-T D^-1 L^-*T P b) = A^-1 b 610*bf2c3715SXin Li dst = m_transpositions.transpose() * dst; 611*bf2c3715SXin Li } 612*bf2c3715SXin Li #endif 613*bf2c3715SXin Li 614*bf2c3715SXin Li /** \internal use x = ldlt_object.solve(x); 615*bf2c3715SXin Li * 616*bf2c3715SXin Li * This is the \em in-place version of solve(). 617*bf2c3715SXin Li * 618*bf2c3715SXin Li * \param bAndX represents both the right-hand side matrix b and result x. 619*bf2c3715SXin Li * 620*bf2c3715SXin Li * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD. 621*bf2c3715SXin Li * 622*bf2c3715SXin Li * This version avoids a copy when the right hand side matrix b is not 623*bf2c3715SXin Li * needed anymore. 624*bf2c3715SXin Li * 625*bf2c3715SXin Li * \sa LDLT::solve(), MatrixBase::ldlt() 626*bf2c3715SXin Li */ 627*bf2c3715SXin Li template<typename MatrixType,int _UpLo> 628*bf2c3715SXin Li template<typename Derived> 629*bf2c3715SXin Li bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const 630*bf2c3715SXin Li { 631*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 632*bf2c3715SXin Li eigen_assert(m_matrix.rows() == bAndX.rows()); 633*bf2c3715SXin Li 634*bf2c3715SXin Li bAndX = this->solve(bAndX); 635*bf2c3715SXin Li 636*bf2c3715SXin Li return true; 637*bf2c3715SXin Li } 638*bf2c3715SXin Li 639*bf2c3715SXin Li /** \returns the matrix represented by the decomposition, 640*bf2c3715SXin Li * i.e., it returns the product: P^T L D L^* P. 641*bf2c3715SXin Li * This function is provided for debug purpose. */ 642*bf2c3715SXin Li template<typename MatrixType, int _UpLo> 643*bf2c3715SXin Li MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const 644*bf2c3715SXin Li { 645*bf2c3715SXin Li eigen_assert(m_isInitialized && "LDLT is not initialized."); 646*bf2c3715SXin Li const Index size = m_matrix.rows(); 647*bf2c3715SXin Li MatrixType res(size,size); 648*bf2c3715SXin Li 649*bf2c3715SXin Li // P 650*bf2c3715SXin Li res.setIdentity(); 651*bf2c3715SXin Li res = transpositionsP() * res; 652*bf2c3715SXin Li // L^* P 653*bf2c3715SXin Li res = matrixU() * res; 654*bf2c3715SXin Li // D(L^*P) 655*bf2c3715SXin Li res = vectorD().real().asDiagonal() * res; 656*bf2c3715SXin Li // L(DL^*P) 657*bf2c3715SXin Li res = matrixL() * res; 658*bf2c3715SXin Li // P^T (LDL^*P) 659*bf2c3715SXin Li res = transpositionsP().transpose() * res; 660*bf2c3715SXin Li 661*bf2c3715SXin Li return res; 662*bf2c3715SXin Li } 663*bf2c3715SXin Li 664*bf2c3715SXin Li /** \cholesky_module 665*bf2c3715SXin Li * \returns the Cholesky decomposition with full pivoting without square root of \c *this 666*bf2c3715SXin Li * \sa MatrixBase::ldlt() 667*bf2c3715SXin Li */ 668*bf2c3715SXin Li template<typename MatrixType, unsigned int UpLo> 669*bf2c3715SXin Li inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo> 670*bf2c3715SXin Li SelfAdjointView<MatrixType, UpLo>::ldlt() const 671*bf2c3715SXin Li { 672*bf2c3715SXin Li return LDLT<PlainObject,UpLo>(m_matrix); 673*bf2c3715SXin Li } 674*bf2c3715SXin Li 675*bf2c3715SXin Li /** \cholesky_module 676*bf2c3715SXin Li * \returns the Cholesky decomposition with full pivoting without square root of \c *this 677*bf2c3715SXin Li * \sa SelfAdjointView::ldlt() 678*bf2c3715SXin Li */ 679*bf2c3715SXin Li template<typename Derived> 680*bf2c3715SXin Li inline const LDLT<typename MatrixBase<Derived>::PlainObject> 681*bf2c3715SXin Li MatrixBase<Derived>::ldlt() const 682*bf2c3715SXin Li { 683*bf2c3715SXin Li return LDLT<PlainObject>(derived()); 684*bf2c3715SXin Li } 685*bf2c3715SXin Li 686*bf2c3715SXin Li } // end namespace Eigen 687*bf2c3715SXin Li 688*bf2c3715SXin Li #endif // EIGEN_LDLT_H 689