LDLT.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
6 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7 // Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
8 //
9 // This Source Code Form is subject to the terms of the Mozilla
10 // Public License v. 2.0. If a copy of the MPL was not distributed
11 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12 
13 #ifndef EIGEN_LDLT_H
14 #define EIGEN_LDLT_H
15 
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21  template<typename MatrixType_, int UpLo_> struct traits<LDLT<MatrixType_, UpLo_> >
22  : traits<MatrixType_>
23  {
24  typedef MatrixXpr XprKind;
25  typedef SolverStorage StorageKind;
26  typedef int StorageIndex;
27  enum { Flags = 0 };
28  };
29 
30  template<typename MatrixType, int UpLo> struct LDLT_Traits;
31 
32  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
34 }
35 
61 template<typename MatrixType_, int UpLo_> class LDLT
62  : public SolverBase<LDLT<MatrixType_, UpLo_> >
63 {
64  public:
65  typedef MatrixType_ MatrixType;
67  friend class SolverBase<LDLT>;
68 
70  enum {
71  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
73  UpLo = UpLo_
74  };
76 
79 
80  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
81 
87  LDLT()
88  : m_matrix(),
91  m_isInitialized(false)
92  {}
93 
100  explicit LDLT(Index size)
101  : m_matrix(size, size),
103  m_temporary(size),
105  m_isInitialized(false)
106  {}
107 
114  template<typename InputType>
115  explicit LDLT(const EigenBase<InputType>& matrix)
116  : m_matrix(matrix.rows(), matrix.cols()),
117  m_transpositions(matrix.rows()),
118  m_temporary(matrix.rows()),
120  m_isInitialized(false)
121  {
122  compute(matrix.derived());
123  }
124 
131  template<typename InputType>
132  explicit LDLT(EigenBase<InputType>& matrix)
133  : m_matrix(matrix.derived()),
134  m_transpositions(matrix.rows()),
135  m_temporary(matrix.rows()),
137  m_isInitialized(false)
138  {
139  compute(matrix.derived());
140  }
141 
145  void setZero()
146  {
147  m_isInitialized = false;
148  }
149 
151  inline typename Traits::MatrixU matrixU() const
152  {
153  eigen_assert(m_isInitialized && "LDLT is not initialized.");
154  return Traits::getU(m_matrix);
155  }
156 
158  inline typename Traits::MatrixL matrixL() const
159  {
160  eigen_assert(m_isInitialized && "LDLT is not initialized.");
161  return Traits::getL(m_matrix);
162  }
163 
166  inline const TranspositionType& transpositionsP() const
167  {
168  eigen_assert(m_isInitialized && "LDLT is not initialized.");
169  return m_transpositions;
170  }
171 
174  {
175  eigen_assert(m_isInitialized && "LDLT is not initialized.");
176  return m_matrix.diagonal();
177  }
178 
180  inline bool isPositive() const
181  {
182  eigen_assert(m_isInitialized && "LDLT is not initialized.");
184  }
185 
187  inline bool isNegative(void) const
188  {
189  eigen_assert(m_isInitialized && "LDLT is not initialized.");
191  }
192 
193  #ifdef EIGEN_PARSED_BY_DOXYGEN
209  template<typename Rhs>
210  inline const Solve<LDLT, Rhs>
211  solve(const MatrixBase<Rhs>& b) const;
212  #endif
213 
214  template<typename Derived>
215  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
216 
217  template<typename InputType>
219 
223  RealScalar rcond() const
224  {
225  eigen_assert(m_isInitialized && "LDLT is not initialized.");
227  }
228 
229  template <typename Derived>
230  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
231 
236  inline const MatrixType& matrixLDLT() const
237  {
238  eigen_assert(m_isInitialized && "LDLT is not initialized.");
239  return m_matrix;
240  }
241 
243 
249  const LDLT& adjoint() const { return *this; }
250 
253 
260  {
261  eigen_assert(m_isInitialized && "LDLT is not initialized.");
262  return m_info;
263  }
264 
265  #ifndef EIGEN_PARSED_BY_DOXYGEN
266  template<typename RhsType, typename DstType>
267  void _solve_impl(const RhsType &rhs, DstType &dst) const;
268 
269  template<bool Conjugate, typename RhsType, typename DstType>
270  void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
271  #endif
272 
273  protected:
274 
276 
277 
284  RealScalar m_l1_norm;
290 };
291 
292 namespace internal {
293 
294 template<int UpLo> struct ldlt_inplace;
295 
296 template<> struct ldlt_inplace<Lower>
297 {
298  template<typename MatrixType, typename TranspositionType, typename Workspace>
299  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
300  {
301  using std::abs;
302  typedef typename MatrixType::Scalar Scalar;
303  typedef typename MatrixType::RealScalar RealScalar;
304  typedef typename TranspositionType::StorageIndex IndexType;
305  eigen_assert(mat.rows()==mat.cols());
306  const Index size = mat.rows();
307  bool found_zero_pivot = false;
308  bool ret = true;
309 
310  if (size <= 1)
311  {
312  transpositions.setIdentity();
313  if(size==0) sign = ZeroSign;
314  else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
315  else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
316  else sign = ZeroSign;
317  return true;
318  }
319 
320  for (Index k = 0; k < size; ++k)
321  {
322  // Find largest diagonal element
323  Index index_of_biggest_in_corner;
324  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
325  index_of_biggest_in_corner += k;
326 
327  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
328  if(k != index_of_biggest_in_corner)
329  {
330  // apply the transposition while taking care to consider only
331  // the lower triangular part
332  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
333  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
334  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
335  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
336  for(Index i=k+1;i<index_of_biggest_in_corner;++i)
337  {
338  Scalar tmp = mat.coeffRef(i,k);
339  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
340  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
341  }
343  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
344  }
345 
346  // partition the matrix:
347  // A00 | - | -
348  // lu = A10 | A11 | -
349  // A20 | A21 | A22
350  Index rs = size - k - 1;
351  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
352  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
353  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
354 
355  if(k>0)
356  {
357  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
358  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
359  if(rs>0)
360  A21.noalias() -= A20 * temp.head(k);
361  }
362 
363  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
364  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
365  // we should only make sure that we do not introduce INF or NaN values.
366  // Remark that LAPACK also uses 0 as the cutoff value.
367  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
368  bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
369 
370  if(k==0 && !pivot_is_valid)
371  {
372  // The entire diagonal is zero, there is nothing more to do
373  // except filling the transpositions, and checking whether the matrix is zero.
374  sign = ZeroSign;
375  for(Index j = 0; j<size; ++j)
376  {
377  transpositions.coeffRef(j) = IndexType(j);
378  ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
379  }
380  return ret;
381  }
382 
383  if((rs>0) && pivot_is_valid)
384  A21 /= realAkk;
385  else if(rs>0)
386  ret = ret && (A21.array()==Scalar(0)).all();
387 
388  if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
389  else if(!pivot_is_valid) found_zero_pivot = true;
390 
391  if (sign == PositiveSemiDef) {
392  if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
393  } else if (sign == NegativeSemiDef) {
394  if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
395  } else if (sign == ZeroSign) {
396  if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
397  else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
398  }
399  }
400 
401  return ret;
402  }
403 
404  // Reference for the algorithm: Davis and Hager, "Multiple Rank
405  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
406  // Trivial rearrangements of their computations (Timothy E. Holy)
407  // allow their algorithm to work for rank-1 updates even if the
408  // original matrix is not of full rank.
409  // Here only rank-1 updates are implemented, to reduce the
410  // requirement for intermediate storage and improve accuracy
411  template<typename MatrixType, typename WDerived>
412  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
413  {
414  using numext::isfinite;
415  typedef typename MatrixType::Scalar Scalar;
416  typedef typename MatrixType::RealScalar RealScalar;
417 
418  const Index size = mat.rows();
419  eigen_assert(mat.cols() == size && w.size()==size);
420 
421  RealScalar alpha = 1;
422 
423  // Apply the update
424  for (Index j = 0; j < size; j++)
425  {
426  // Check for termination due to an original decomposition of low-rank
427  if (!(isfinite)(alpha))
428  break;
429 
430  // Update the diagonal terms
431  RealScalar dj = numext::real(mat.coeff(j,j));
432  Scalar wj = w.coeff(j);
433  RealScalar swj2 = sigma*numext::abs2(wj);
434  RealScalar gamma = dj*alpha + swj2;
435 
436  mat.coeffRef(j,j) += swj2/alpha;
437  alpha += swj2/dj;
438 
439 
440  // Update the terms of L
441  Index rs = size-j-1;
442  w.tail(rs) -= wj * mat.col(j).tail(rs);
443  if(!numext::is_exactly_zero(gamma))
444  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
445  }
446  return true;
447  }
448 
449  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
450  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
451  {
452  // Apply the permutation to the input w
453  tmp = transpositions * w;
454 
455  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
456  }
457 };
458 
459 template<> struct ldlt_inplace<Upper>
460 {
461  template<typename MatrixType, typename TranspositionType, typename Workspace>
462  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
463  {
464  Transpose<MatrixType> matt(mat);
465  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
466  }
467 
468  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
469  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
470  {
471  Transpose<MatrixType> matt(mat);
472  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
473  }
474 };
475 
476 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
477 {
478  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
479  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
480  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
481  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
482 };
483 
484 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
485 {
486  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
487  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
488  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
489  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
490 };
491 
492 } // end namespace internal
493 
496 template<typename MatrixType, int UpLo_>
497 template<typename InputType>
499 {
500  eigen_assert(a.rows()==a.cols());
501  const Index size = a.rows();
502 
503  m_matrix = a.derived();
504 
505  // Compute matrix L1 norm = max abs column sum.
506  m_l1_norm = RealScalar(0);
507  // TODO move this code to SelfAdjointView
508  for (Index col = 0; col < size; ++col) {
509  RealScalar abs_col_sum;
510  if (UpLo_ == Lower)
511  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
512  else
513  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
514  if (abs_col_sum > m_l1_norm)
515  m_l1_norm = abs_col_sum;
516  }
517 
519  m_isInitialized = false;
522 
523  m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
524 
525  m_isInitialized = true;
526  return *this;
527 }
528 
534 template<typename MatrixType, int UpLo_>
535 template<typename Derived>
537 {
538  typedef typename TranspositionType::StorageIndex IndexType;
539  const Index size = w.rows();
540  if (m_isInitialized)
541  {
542  eigen_assert(m_matrix.rows()==size);
543  }
544  else
545  {
546  m_matrix.resize(size,size);
547  m_matrix.setZero();
549  for (Index i = 0; i < size; i++)
550  m_transpositions.coeffRef(i) = IndexType(i);
553  m_isInitialized = true;
554  }
555 
556  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
557 
558  return *this;
559 }
560 
561 #ifndef EIGEN_PARSED_BY_DOXYGEN
562 template<typename MatrixType_, int UpLo_>
563 template<typename RhsType, typename DstType>
564 void LDLT<MatrixType_,UpLo_>::_solve_impl(const RhsType &rhs, DstType &dst) const
565 {
566  _solve_impl_transposed<true>(rhs, dst);
567 }
568 
569 template<typename MatrixType_,int UpLo_>
570 template<bool Conjugate, typename RhsType, typename DstType>
571 void LDLT<MatrixType_,UpLo_>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
572 {
573  // dst = P b
574  dst = m_transpositions * rhs;
575 
576  // dst = L^-1 (P b)
577  // dst = L^-*T (P b)
578  matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
579 
580  // dst = D^-* (L^-1 P b)
581  // dst = D^-1 (L^-*T P b)
582  // more precisely, use pseudo-inverse of D (see bug 241)
583  using std::abs;
584  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
585  // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
586  // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
587  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
588  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
589  // diagonal element is not well justified and leads to numerical issues in some cases.
590  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
591  // Using numeric_limits::min() gives us more robustness to denormals.
592  RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
593  for (Index i = 0; i < vecD.size(); ++i)
594  {
595  if(abs(vecD(i)) > tolerance)
596  dst.row(i) /= vecD(i);
597  else
598  dst.row(i).setZero();
599  }
600 
601  // dst = L^-* (D^-* L^-1 P b)
602  // dst = L^-T (D^-1 L^-*T P b)
603  matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
604 
605  // dst = P^T (L^-* D^-* L^-1 P b) = A^-1 b
606  // dst = P^-T (L^-T D^-1 L^-*T P b) = A^-1 b
607  dst = m_transpositions.transpose() * dst;
608 }
609 #endif
610 
624 template<typename MatrixType,int UpLo_>
625 template<typename Derived>
627 {
628  eigen_assert(m_isInitialized && "LDLT is not initialized.");
629  eigen_assert(m_matrix.rows() == bAndX.rows());
630 
631  bAndX = this->solve(bAndX);
632 
633  return true;
634 }
635 
639 template<typename MatrixType, int UpLo_>
641 {
642  eigen_assert(m_isInitialized && "LDLT is not initialized.");
643  const Index size = m_matrix.rows();
645 
646  // P
647  res.setIdentity();
648  res = transpositionsP() * res;
649  // L^* P
650  res = matrixU() * res;
651  // D(L^*P)
652  res = vectorD().real().asDiagonal() * res;
653  // L(DL^*P)
654  res = matrixL() * res;
655  // P^T (LDL^*P)
657 
658  return res;
659 }
660 
665 template<typename MatrixType, unsigned int UpLo>
668 {
670 }
671 
676 template<typename Derived>
679 {
680  return LDLT<PlainObject>(derived());
681 }
682 
683 } // end namespace Eigen
684 
685 #endif // EIGEN_LDLT_H
Matrix3f m
const AbsReturnType abs() const
Array< int, 3, 1 > b
ColXpr col(Index i)
This is the const version of col().
RealReturnType real() const
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
Definition: Macros.h:1149
#define EIGEN_NOEXCEPT
Definition: Macros.h:1260
#define EIGEN_CONSTEXPR
Definition: Macros.h:747
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
RowVector3d w
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:81
Matrix< float, 1, Dynamic > MatrixType
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:107
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:67
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:63
MatrixType_ MatrixType
Definition: LDLT.h:65
LDLT & compute(const EigenBase< InputType > &matrix)
const TranspositionType & transpositionsP() const
Definition: LDLT.h:166
@ MaxColsAtCompileTime
Definition: LDLT.h:72
@ MaxRowsAtCompileTime
Definition: LDLT.h:71
@ UpLo
Definition: LDLT.h:73
internal::LDLT_Traits< MatrixType, UpLo > Traits
Definition: LDLT.h:80
Traits::MatrixL matrixL() const
Definition: LDLT.h:158
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:173
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
ComputationInfo m_info
Definition: LDLT.h:289
internal::SignMatrix m_sign
Definition: LDLT.h:287
const MatrixType & matrixLDLT() const
Definition: LDLT.h:236
bool m_isInitialized
Definition: LDLT.h:288
void setZero()
Definition: LDLT.h:145
TranspositionType m_transpositions
Definition: LDLT.h:285
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
Definition: LDLT.h:77
RealScalar rcond() const
Definition: LDLT.h:223
MatrixType m_matrix
Definition: LDLT.h:283
bool isNegative(void) const
Definition: LDLT.h:187
LDLT()
Default Constructor.
Definition: LDLT.h:87
bool isPositive() const
Definition: LDLT.h:180
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: LDLT.h:251
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:115
SolverBase< LDLT > Base
Definition: LDLT.h:66
RealScalar m_l1_norm
Definition: LDLT.h:284
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:100
LDLT & rankUpdate(const MatrixBase< Derived > &w, const RealScalar &alpha=1)
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:132
Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
Definition: LDLT.h:75
TmpMatrixType m_temporary
Definition: LDLT.h:286
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: LDLT.h:252
MatrixType reconstructedMatrix() const
Definition: LDLT.h:640
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:259
Traits::MatrixU matrixU() const
Definition: LDLT.h:151
bool solveInPlace(MatrixBase< Derived > &bAndX) const
Definition: LDLT.h:626
const LDLT & adjoint() const
Definition: LDLT.h:249
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationType
Definition: LDLT.h:78
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:678
constexpr void resize(Index rows, Index cols)
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:667
Pseudo expression representing a solving operation.
Definition: Solve.h:65
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:71
internal::traits< LDLT< MatrixType_, UpLo_ > >::Scalar Scalar
Definition: SolverBase.h:75
LDLT< MatrixType_, UpLo_ > & derived()
Definition: EigenBase.h:48
Transpose< TranspositionsBase > transpose() const
void resize(Index newSize)
StorageIndex & coeffRef(Index i)
static const Eigen::internal::all_t all
ComputationInfo
Definition: Constants.h:444
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
@ PositiveSemiDef
Definition: LDLT.h:33
@ NegativeSemiDef
Definition: LDLT.h:33
@ Indefinite
Definition: LDLT.h:33
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:788
bool abs2(bool x)
EIGEN_ALWAYS_INLINE bool() isfinite(const Eigen::bfloat16 &h)
Definition: BFloat16.h:790
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Derived & derived()
Definition: EigenBase.h:48
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231
std::ptrdiff_t j