LDLT.h 23.9 KB
Newer Older
LM's avatar
LM committed
1 2 3
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
Don Gagne's avatar
Don Gagne committed
4
// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
LM's avatar
LM committed
5 6
// Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
Don Gagne's avatar
Don Gagne committed
7
// Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
LM's avatar
LM committed
8
//
Don Gagne's avatar
Don Gagne committed
9 10 11
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
12 13 14 15

#ifndef EIGEN_LDLT_H
#define EIGEN_LDLT_H

16
namespace Eigen {
Don Gagne's avatar
Don Gagne committed
17

LM's avatar
LM committed
18
namespace internal {
19 20 21 22
  template<typename MatrixType, int UpLo> struct LDLT_Traits;

  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
  enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
LM's avatar
LM committed
23 24
}

Don Gagne's avatar
Don Gagne committed
25
/** \ingroup Cholesky_Module
LM's avatar
LM committed
26 27 28 29 30
  *
  * \class LDLT
  *
  * \brief Robust Cholesky decomposition of a matrix with pivoting
  *
31 32
  * \tparam _MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
  * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
Don Gagne's avatar
Don Gagne committed
33
  *             The other triangular part won't be read.
LM's avatar
LM committed
34 35 36 37 38 39 40 41 42 43
  *
  * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
  * matrix \f$ A \f$ such that \f$ A =  P^TLDL^*P \f$, where P is a permutation matrix, L
  * is lower triangular with a unit diagonal and D is a diagonal matrix.
  *
  * The decomposition uses pivoting to ensure stability, so that L will have
  * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
  * on D also stabilizes the computation.
  *
  * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
Don Gagne's avatar
Don Gagne committed
44
  * decomposition to determine whether a system of equations has a solution.
LM's avatar
LM committed
45
  *
46 47 48
  * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
  * 
  * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT
LM's avatar
LM committed
49 50 51 52 53 54 55 56 57 58 59 60 61 62
  */
template<typename _MatrixType, int _UpLo> class LDLT
{
  public:
    typedef _MatrixType MatrixType;
    enum {
      RowsAtCompileTime = MatrixType::RowsAtCompileTime,
      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
      MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
      UpLo = _UpLo
    };
    typedef typename MatrixType::Scalar Scalar;
    typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
63 64 65
    typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
    typedef typename MatrixType::StorageIndex StorageIndex;
    typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> TmpMatrixType;
LM's avatar
LM committed
66 67 68 69 70 71 72 73 74 75 76

    typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
    typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;

    typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;

    /** \brief Default Constructor.
      *
      * The default constructor is useful in cases in which the user intends to
      * perform decompositions via LDLT::compute(const MatrixType&).
      */
77 78 79
    LDLT()
      : m_matrix(),
        m_transpositions(),
80
        m_sign(internal::ZeroSign),
81
        m_isInitialized(false)
82
    {}
LM's avatar
LM committed
83 84 85 86 87 88 89

    /** \brief Default Constructor with memory preallocation
      *
      * Like the default constructor but with preallocation of the internal data
      * according to the specified problem \a size.
      * \sa LDLT()
      */
90
    explicit LDLT(Index size)
LM's avatar
LM committed
91 92 93
      : m_matrix(size, size),
        m_transpositions(size),
        m_temporary(size),
94
        m_sign(internal::ZeroSign),
LM's avatar
LM committed
95 96 97
        m_isInitialized(false)
    {}

Don Gagne's avatar
Don Gagne committed
98 99 100
    /** \brief Constructor with decomposition
      *
      * This calculates the decomposition for the input \a matrix.
101
      *
Don Gagne's avatar
Don Gagne committed
102 103
      * \sa LDLT(Index size)
      */
104 105
    template<typename InputType>
    explicit LDLT(const EigenBase<InputType>& matrix)
LM's avatar
LM committed
106 107 108
      : m_matrix(matrix.rows(), matrix.cols()),
        m_transpositions(matrix.rows()),
        m_temporary(matrix.rows()),
109
        m_sign(internal::ZeroSign),
LM's avatar
LM committed
110 111
        m_isInitialized(false)
    {
112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129
      compute(matrix.derived());
    }

    /** \brief Constructs a LDLT factorization from a given matrix
      *
      * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref.
      *
      * \sa LDLT(const EigenBase&)
      */
    template<typename InputType>
    explicit LDLT(EigenBase<InputType>& matrix)
      : m_matrix(matrix.derived()),
        m_transpositions(matrix.rows()),
        m_temporary(matrix.rows()),
        m_sign(internal::ZeroSign),
        m_isInitialized(false)
    {
      compute(matrix.derived());
LM's avatar
LM committed
130 131
    }

Don Gagne's avatar
Don Gagne committed
132 133 134 135 136 137 138 139
    /** Clear any existing decomposition
     * \sa rankUpdate(w,sigma)
     */
    void setZero()
    {
      m_isInitialized = false;
    }

LM's avatar
LM committed
140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162
    /** \returns a view of the upper triangular matrix U */
    inline typename Traits::MatrixU matrixU() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return Traits::getU(m_matrix);
    }

    /** \returns a view of the lower triangular matrix L */
    inline typename Traits::MatrixL matrixL() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return Traits::getL(m_matrix);
    }

    /** \returns the permutation matrix P as a transposition sequence.
      */
    inline const TranspositionType& transpositionsP() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return m_transpositions;
    }

    /** \returns the coefficients of the diagonal matrix D */
Don Gagne's avatar
Don Gagne committed
163
    inline Diagonal<const MatrixType> vectorD() const
LM's avatar
LM committed
164 165 166 167 168 169
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return m_matrix.diagonal();
    }

    /** \returns true if the matrix is positive (semidefinite) */
Don Gagne's avatar
Don Gagne committed
170
    inline bool isPositive() const
LM's avatar
LM committed
171 172
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
173
      return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
LM's avatar
LM committed
174 175 176 177 178 179
    }

    /** \returns true if the matrix is negative (semidefinite) */
    inline bool isNegative(void) const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
180
      return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
LM's avatar
LM committed
181 182 183
    }

    /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
Don Gagne's avatar
Don Gagne committed
184 185
      *
      * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
LM's avatar
LM committed
186 187 188
      *
      * \note_about_checking_solutions
      *
Don Gagne's avatar
Don Gagne committed
189
      * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
190
      * 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$,
Don Gagne's avatar
Don Gagne committed
191 192 193 194 195
      * \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
      * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
      * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
      * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
      *
196
      * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt()
LM's avatar
LM committed
197 198
      */
    template<typename Rhs>
199
    inline const Solve<LDLT, Rhs>
LM's avatar
LM committed
200 201 202 203 204
    solve(const MatrixBase<Rhs>& b) const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      eigen_assert(m_matrix.rows()==b.rows()
                && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
205
      return Solve<LDLT, Rhs>(*this, b.derived());
LM's avatar
LM committed
206 207 208 209 210
    }

    template<typename Derived>
    bool solveInPlace(MatrixBase<Derived> &bAndX) const;

211 212 213 214 215 216 217 218 219 220 221
    template<typename InputType>
    LDLT& compute(const EigenBase<InputType>& matrix);

    /** \returns an estimate of the reciprocal condition number of the matrix of
     *  which \c *this is the LDLT decomposition.
     */
    RealScalar rcond() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return internal::rcond_estimate_helper(m_l1_norm, *this);
    }
LM's avatar
LM committed
222

Don Gagne's avatar
Don Gagne committed
223 224 225
    template <typename Derived>
    LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);

LM's avatar
LM committed
226 227 228 229 230 231 232 233 234 235 236 237
    /** \returns the internal LDLT decomposition matrix
      *
      * TODO: document the storage layout
      */
    inline const MatrixType& matrixLDLT() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
      return m_matrix;
    }

    MatrixType reconstructedMatrix() const;

238 239 240 241 242 243 244
    /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.
      *
      * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:
      * \code x = decomposition.adjoint().solve(b) \endcode
      */
    const LDLT& adjoint() const { return *this; };

LM's avatar
LM committed
245 246 247
    inline Index rows() const { return m_matrix.rows(); }
    inline Index cols() const { return m_matrix.cols(); }

Don Gagne's avatar
Don Gagne committed
248 249 250
    /** \brief Reports whether previous computation was successful.
      *
      * \returns \c Success if computation was succesful,
251
      *          \c NumericalIssue if the factorization failed because of a zero pivot.
Don Gagne's avatar
Don Gagne committed
252 253 254 255
      */
    ComputationInfo info() const
    {
      eigen_assert(m_isInitialized && "LDLT is not initialized.");
256
      return m_info;
Don Gagne's avatar
Don Gagne committed
257 258
    }

259 260 261 262 263 264
    #ifndef EIGEN_PARSED_BY_DOXYGEN
    template<typename RhsType, typename DstType>
    EIGEN_DEVICE_FUNC
    void _solve_impl(const RhsType &rhs, DstType &dst) const;
    #endif

LM's avatar
LM committed
265
  protected:
266

267 268 269 270
    static void check_template_parameters()
    {
      EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
    }
LM's avatar
LM committed
271 272 273 274 275 276 277 278

    /** \internal
      * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
      * The strict upper part is used during the decomposition, the strict lower
      * part correspond to the coefficients of L (its diagonal is equal to 1 and
      * is not stored), and the diagonal entries correspond to D.
      */
    MatrixType m_matrix;
279
    RealScalar m_l1_norm;
LM's avatar
LM committed
280 281
    TranspositionType m_transpositions;
    TmpMatrixType m_temporary;
282
    internal::SignMatrix m_sign;
LM's avatar
LM committed
283
    bool m_isInitialized;
284
    ComputationInfo m_info;
LM's avatar
LM committed
285 286 287 288 289 290 291 292 293
};

namespace internal {

template<int UpLo> struct ldlt_inplace;

template<> struct ldlt_inplace<Lower>
{
  template<typename MatrixType, typename TranspositionType, typename Workspace>
294
  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
LM's avatar
LM committed
295
  {
Don Gagne's avatar
Don Gagne committed
296
    using std::abs;
LM's avatar
LM committed
297 298
    typedef typename MatrixType::Scalar Scalar;
    typedef typename MatrixType::RealScalar RealScalar;
299
    typedef typename TranspositionType::StorageIndex IndexType;
LM's avatar
LM committed
300 301
    eigen_assert(mat.rows()==mat.cols());
    const Index size = mat.rows();
302 303
    bool found_zero_pivot = false;
    bool ret = true;
LM's avatar
LM committed
304 305 306 307

    if (size <= 1)
    {
      transpositions.setIdentity();
308 309 310
      if(size==0) sign = ZeroSign;
      else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
      else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
311
      else sign = ZeroSign;
LM's avatar
LM committed
312 313 314 315 316 317 318
      return true;
    }

    for (Index k = 0; k < size; ++k)
    {
      // Find largest diagonal element
      Index index_of_biggest_in_corner;
319
      mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
LM's avatar
LM committed
320 321
      index_of_biggest_in_corner += k;

322
      transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
LM's avatar
LM committed
323 324 325 326 327 328 329 330
      if(k != index_of_biggest_in_corner)
      {
        // apply the transposition while taking care to consider only
        // the lower triangular part
        Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
        mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
        mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
        std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
331
        for(Index i=k+1;i<index_of_biggest_in_corner;++i)
LM's avatar
LM committed
332 333
        {
          Scalar tmp = mat.coeffRef(i,k);
Don Gagne's avatar
Don Gagne committed
334 335
          mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
          mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
LM's avatar
LM committed
336 337
        }
        if(NumTraits<Scalar>::IsComplex)
Don Gagne's avatar
Don Gagne committed
338
          mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
LM's avatar
LM committed
339 340 341 342 343 344 345 346 347 348 349 350 351
      }

      // partition the matrix:
      //       A00 |  -  |  -
      // lu  = A10 | A11 |  -
      //       A20 | A21 | A22
      Index rs = size - k - 1;
      Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
      Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
      Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);

      if(k>0)
      {
352
        temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
LM's avatar
LM committed
353 354 355 356
        mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
        if(rs>0)
          A21.noalias() -= A20 * temp.head(k);
      }
357

358
      // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
359 360 361
      // was smaller than the cutoff value. However, since LDLT is not rank-revealing
      // we should only make sure that we do not introduce INF or NaN values.
      // Remark that LAPACK also uses 0 as the cutoff value.
362
      RealScalar realAkk = numext::real(mat.coeffRef(k,k));
363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378
      bool pivot_is_valid = (abs(realAkk) > RealScalar(0));

      if(k==0 && !pivot_is_valid)
      {
        // The entire diagonal is zero, there is nothing more to do
        // except filling the transpositions, and checking whether the matrix is zero.
        sign = ZeroSign;
        for(Index j = 0; j<size; ++j)
        {
          transpositions.coeffRef(j) = IndexType(j);
          ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
        }
        return ret;
      }

      if((rs>0) && pivot_is_valid)
379
        A21 /= realAkk;
380 381 382 383 384
      else if(rs>0)
        ret = ret && (A21.array()==Scalar(0)).all();

      if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
      else if(!pivot_is_valid) found_zero_pivot = true;
385 386

      if (sign == PositiveSemiDef) {
387
        if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
388
      } else if (sign == NegativeSemiDef) {
389
        if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
390
      } else if (sign == ZeroSign) {
391 392
        if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
        else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
Don Gagne's avatar
Don Gagne committed
393
      }
LM's avatar
LM committed
394 395
    }

396
    return ret;
LM's avatar
LM committed
397
  }
Don Gagne's avatar
Don Gagne committed
398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451

  // Reference for the algorithm: Davis and Hager, "Multiple Rank
  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
  // Trivial rearrangements of their computations (Timothy E. Holy)
  // allow their algorithm to work for rank-1 updates even if the
  // original matrix is not of full rank.
  // Here only rank-1 updates are implemented, to reduce the
  // requirement for intermediate storage and improve accuracy
  template<typename MatrixType, typename WDerived>
  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
  {
    using numext::isfinite;
    typedef typename MatrixType::Scalar Scalar;
    typedef typename MatrixType::RealScalar RealScalar;

    const Index size = mat.rows();
    eigen_assert(mat.cols() == size && w.size()==size);

    RealScalar alpha = 1;

    // Apply the update
    for (Index j = 0; j < size; j++)
    {
      // Check for termination due to an original decomposition of low-rank
      if (!(isfinite)(alpha))
        break;

      // Update the diagonal terms
      RealScalar dj = numext::real(mat.coeff(j,j));
      Scalar wj = w.coeff(j);
      RealScalar swj2 = sigma*numext::abs2(wj);
      RealScalar gamma = dj*alpha + swj2;

      mat.coeffRef(j,j) += swj2/alpha;
      alpha += swj2/dj;


      // Update the terms of L
      Index rs = size-j-1;
      w.tail(rs) -= wj * mat.col(j).tail(rs);
      if(gamma != 0)
        mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
    }
    return true;
  }

  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
  {
    // Apply the permutation to the input w
    tmp = transpositions * w;

    return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
  }
LM's avatar
LM committed
452 453 454 455 456
};

template<> struct ldlt_inplace<Upper>
{
  template<typename MatrixType, typename TranspositionType, typename Workspace>
457
  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
LM's avatar
LM committed
458 459 460 461
  {
    Transpose<MatrixType> matt(mat);
    return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
  }
Don Gagne's avatar
Don Gagne committed
462 463 464 465 466 467 468

  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
  {
    Transpose<MatrixType> matt(mat);
    return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
  }
LM's avatar
LM committed
469 470 471 472
};

template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
{
Don Gagne's avatar
Don Gagne committed
473 474
  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
475 476
  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
LM's avatar
LM committed
477 478 479 480
};

template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
{
Don Gagne's avatar
Don Gagne committed
481 482
  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
483 484
  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
LM's avatar
LM committed
485 486 487 488 489 490 491
};

} // end namespace internal

/** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
  */
template<typename MatrixType, int _UpLo>
492 493
template<typename InputType>
LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
LM's avatar
LM committed
494
{
495
  check_template_parameters();
496

LM's avatar
LM committed
497 498 499
  eigen_assert(a.rows()==a.cols());
  const Index size = a.rows();

500 501 502 503 504 505 506 507 508 509 510 511 512 513
  m_matrix = a.derived();

  // Compute matrix L1 norm = max abs column sum.
  m_l1_norm = RealScalar(0);
  // TODO move this code to SelfAdjointView
  for (Index col = 0; col < size; ++col) {
    RealScalar abs_col_sum;
    if (_UpLo == Lower)
      abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
    else
      abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
    if (abs_col_sum > m_l1_norm)
      m_l1_norm = abs_col_sum;
  }
LM's avatar
LM committed
514 515 516 517

  m_transpositions.resize(size);
  m_isInitialized = false;
  m_temporary.resize(size);
518
  m_sign = internal::ZeroSign;
LM's avatar
LM committed
519

520
  m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
LM's avatar
LM committed
521 522 523 524 525

  m_isInitialized = true;
  return *this;
}

Don Gagne's avatar
Don Gagne committed
526 527 528 529 530 531 532
/** Update the LDLT decomposition:  given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
 * \param w a vector to be incorporated into the decomposition.
 * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
 * \sa setZero()
  */
template<typename MatrixType, int _UpLo>
template<typename Derived>
533
LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
Don Gagne's avatar
Don Gagne committed
534
{
535
  typedef typename TranspositionType::StorageIndex IndexType;
Don Gagne's avatar
Don Gagne committed
536 537 538 539 540 541
  const Index size = w.rows();
  if (m_isInitialized)
  {
    eigen_assert(m_matrix.rows()==size);
  }
  else
542
  {
Don Gagne's avatar
Don Gagne committed
543 544 545 546
    m_matrix.resize(size,size);
    m_matrix.setZero();
    m_transpositions.resize(size);
    for (Index i = 0; i < size; i++)
547
      m_transpositions.coeffRef(i) = IndexType(i);
Don Gagne's avatar
Don Gagne committed
548
    m_temporary.resize(size);
549
    m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
Don Gagne's avatar
Don Gagne committed
550 551 552 553 554 555 556 557
    m_isInitialized = true;
  }

  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);

  return *this;
}

558 559 560 561
#ifndef EIGEN_PARSED_BY_DOXYGEN
template<typename _MatrixType, int _UpLo>
template<typename RhsType, typename DstType>
void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
LM's avatar
LM committed
562
{
563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583
  eigen_assert(rhs.rows() == rows());
  // dst = P b
  dst = m_transpositions * rhs;

  // dst = L^-1 (P b)
  matrixL().solveInPlace(dst);

  // dst = D^-1 (L^-1 P b)
  // more precisely, use pseudo-inverse of D (see bug 241)
  using std::abs;
  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
  // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
  // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
  // diagonal element is not well justified and leads to numerical issues in some cases.
  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
  // Using numeric_limits::min() gives us more robustness to denormals.
  RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();

  for (Index i = 0; i < vecD.size(); ++i)
LM's avatar
LM committed
584
  {
585 586 587 588 589
    if(abs(vecD(i)) > tolerance)
      dst.row(i) /= vecD(i);
    else
      dst.row(i).setZero();
  }
LM's avatar
LM committed
590

591 592
  // dst = L^-T (D^-1 L^-1 P b)
  matrixU().solveInPlace(dst);
LM's avatar
LM committed
593

594 595
  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
  dst = m_transpositions.transpose() * dst;
LM's avatar
LM committed
596
}
597
#endif
LM's avatar
LM committed
598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616

/** \internal use x = ldlt_object.solve(x);
  *
  * This is the \em in-place version of solve().
  *
  * \param bAndX represents both the right-hand side matrix b and result x.
  *
  * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
  *
  * This version avoids a copy when the right hand side matrix b is not
  * needed anymore.
  *
  * \sa LDLT::solve(), MatrixBase::ldlt()
  */
template<typename MatrixType,int _UpLo>
template<typename Derived>
bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
{
  eigen_assert(m_isInitialized && "LDLT is not initialized.");
Don Gagne's avatar
Don Gagne committed
617
  eigen_assert(m_matrix.rows() == bAndX.rows());
LM's avatar
LM committed
618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639

  bAndX = this->solve(bAndX);

  return true;
}

/** \returns the matrix represented by the decomposition,
 * i.e., it returns the product: P^T L D L^* P.
 * This function is provided for debug purpose. */
template<typename MatrixType, int _UpLo>
MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
{
  eigen_assert(m_isInitialized && "LDLT is not initialized.");
  const Index size = m_matrix.rows();
  MatrixType res(size,size);

  // P
  res.setIdentity();
  res = transpositionsP() * res;
  // L^* P
  res = matrixU() * res;
  // D(L^*P)
640
  res = vectorD().real().asDiagonal() * res;
LM's avatar
LM committed
641 642 643 644 645 646 647 648 649 650
  // L(DL^*P)
  res = matrixL() * res;
  // P^T (LDL^*P)
  res = transpositionsP().transpose() * res;

  return res;
}

/** \cholesky_module
  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
651
  * \sa MatrixBase::ldlt()
LM's avatar
LM committed
652 653 654 655 656 657 658 659 660 661
  */
template<typename MatrixType, unsigned int UpLo>
inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
SelfAdjointView<MatrixType, UpLo>::ldlt() const
{
  return LDLT<PlainObject,UpLo>(m_matrix);
}

/** \cholesky_module
  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
662
  * \sa SelfAdjointView::ldlt()
LM's avatar
LM committed
663 664 665 666 667 668 669 670
  */
template<typename Derived>
inline const LDLT<typename MatrixBase<Derived>::PlainObject>
MatrixBase<Derived>::ldlt() const
{
  return LDLT<PlainObject>(derived());
}

Don Gagne's avatar
Don Gagne committed
671 672
} // end namespace Eigen

LM's avatar
LM committed
673
#endif // EIGEN_LDLT_H