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 namespace Eigen {
17 
18 namespace internal {
19  template<typename MatrixType, int UpLo> struct LDLT_Traits;
20 
21  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
23 }
24 
50 template<typename _MatrixType, int _UpLo> class LDLT
51 {
52  public:
53  typedef _MatrixType MatrixType;
54  enum {
55  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
58  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
59  UpLo = _UpLo
60  };
61  typedef typename MatrixType::Scalar Scalar;
63  typedef Eigen::Index Index;
64  typedef typename MatrixType::StorageIndex StorageIndex;
66 
69 
71 
77  LDLT()
78  : m_matrix(),
79  m_transpositions(),
80  m_sign(internal::ZeroSign),
81  m_isInitialized(false)
82  {}
83 
90  explicit LDLT(Index size)
91  : m_matrix(size, size),
92  m_transpositions(size),
93  m_temporary(size),
94  m_sign(internal::ZeroSign),
95  m_isInitialized(false)
96  {}
97 
104  template<typename InputType>
106  : m_matrix(matrix.rows(), matrix.cols()),
107  m_transpositions(matrix.rows()),
108  m_temporary(matrix.rows()),
109  m_sign(internal::ZeroSign),
110  m_isInitialized(false)
111  {
112  compute(matrix.derived());
113  }
114 
121  template<typename InputType>
123  : m_matrix(matrix.derived()),
124  m_transpositions(matrix.rows()),
125  m_temporary(matrix.rows()),
126  m_sign(internal::ZeroSign),
127  m_isInitialized(false)
128  {
129  compute(matrix.derived());
130  }
131 
135  void setZero()
136  {
137  m_isInitialized = false;
138  }
139 
141  inline typename Traits::MatrixU matrixU() const
142  {
143  eigen_assert(m_isInitialized && "LDLT is not initialized.");
144  return Traits::getU(m_matrix);
145  }
146 
148  inline typename Traits::MatrixL matrixL() const
149  {
150  eigen_assert(m_isInitialized && "LDLT is not initialized.");
151  return Traits::getL(m_matrix);
152  }
153 
156  inline const TranspositionType& transpositionsP() const
157  {
158  eigen_assert(m_isInitialized && "LDLT is not initialized.");
159  return m_transpositions;
160  }
161 
164  {
165  eigen_assert(m_isInitialized && "LDLT is not initialized.");
166  return m_matrix.diagonal();
167  }
168 
170  inline bool isPositive() const
171  {
172  eigen_assert(m_isInitialized && "LDLT is not initialized.");
173  return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
174  }
175 
177  inline bool isNegative(void) const
178  {
179  eigen_assert(m_isInitialized && "LDLT is not initialized.");
180  return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
181  }
182 
198  template<typename Rhs>
199  inline const Solve<LDLT, Rhs>
200  solve(const MatrixBase<Rhs>& b) const
201  {
202  eigen_assert(m_isInitialized && "LDLT is not initialized.");
203  eigen_assert(m_matrix.rows()==b.rows()
204  && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
205  return Solve<LDLT, Rhs>(*this, b.derived());
206  }
207 
208  template<typename Derived>
209  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
210 
211  template<typename InputType>
213 
217  RealScalar rcond() const
218  {
219  eigen_assert(m_isInitialized && "LDLT is not initialized.");
220  return internal::rcond_estimate_helper(m_l1_norm, *this);
221  }
222 
223  template <typename Derived>
224  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
225 
230  inline const MatrixType& matrixLDLT() const
231  {
232  eigen_assert(m_isInitialized && "LDLT is not initialized.");
233  return m_matrix;
234  }
235 
236  MatrixType reconstructedMatrix() const;
237 
243  const LDLT& adjoint() const { return *this; };
244 
245  inline Index rows() const { return m_matrix.rows(); }
246  inline Index cols() const { return m_matrix.cols(); }
247 
254  {
255  eigen_assert(m_isInitialized && "LDLT is not initialized.");
256  return m_info;
257  }
258 
259  #ifndef EIGEN_PARSED_BY_DOXYGEN
260  template<typename RhsType, typename DstType>
261  EIGEN_DEVICE_FUNC
262  void _solve_impl(const RhsType &rhs, DstType &dst) const;
263  #endif
264 
265  protected:
266 
268  {
270  }
271 
278  MatrixType m_matrix;
279  RealScalar m_l1_norm;
280  TranspositionType m_transpositions;
281  TmpMatrixType m_temporary;
285 };
286 
287 namespace internal {
288 
289 template<int UpLo> struct ldlt_inplace;
290 
291 template<> struct ldlt_inplace<Lower>
292 {
293  template<typename MatrixType, typename TranspositionType, typename Workspace>
294  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
295  {
296  using std::abs;
297  typedef typename MatrixType::Scalar Scalar;
298  typedef typename MatrixType::RealScalar RealScalar;
299  typedef typename TranspositionType::StorageIndex IndexType;
300  eigen_assert(mat.rows()==mat.cols());
301  const Index size = mat.rows();
302  bool found_zero_pivot = false;
303  bool ret = true;
304 
305  if (size <= 1)
306  {
307  transpositions.setIdentity();
308  if(size==0) sign = ZeroSign;
309  else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
310  else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
311  else sign = ZeroSign;
312  return true;
313  }
314 
315  for (Index k = 0; k < size; ++k)
316  {
317  // Find largest diagonal element
318  Index index_of_biggest_in_corner;
319  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
320  index_of_biggest_in_corner += k;
321 
322  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
323  if(k != index_of_biggest_in_corner)
324  {
325  // apply the transposition while taking care to consider only
326  // the lower triangular part
327  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
328  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
329  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
330  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)
332  {
333  Scalar tmp = mat.coeffRef(i,k);
334  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
335  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
336  }
338  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
339  }
340 
341  // partition the matrix:
342  // A00 | - | -
343  // lu = A10 | A11 | -
344  // A20 | A21 | A22
345  Index rs = size - k - 1;
346  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
347  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
348  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
349 
350  if(k>0)
351  {
352  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
353  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
354  if(rs>0)
355  A21.noalias() -= A20 * temp.head(k);
356  }
357 
358  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
359  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
360  // we should only make sure that we do not introduce INF or NaN values.
361  // Remark that LAPACK also uses 0 as the cutoff value.
362  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
363  bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
364 
365  if(k==0 && !pivot_is_valid)
366  {
367  // The entire diagonal is zero, there is nothing more to do
368  // except filling the transpositions, and checking whether the matrix is zero.
369  sign = ZeroSign;
370  for(Index j = 0; j<size; ++j)
371  {
372  transpositions.coeffRef(j) = IndexType(j);
373  ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
374  }
375  return ret;
376  }
377 
378  if((rs>0) && pivot_is_valid)
379  A21 /= realAkk;
380  else if(rs>0)
381  ret = ret && (A21.array()==Scalar(0)).all();
382 
383  if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
384  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  if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
392  else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
393  }
394  }
395 
396  return ret;
397  }
398 
399  // Reference for the algorithm: Davis and Hager, "Multiple Rank
400  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
401  // Trivial rearrangements of their computations (Timothy E. Holy)
402  // allow their algorithm to work for rank-1 updates even if the
403  // original matrix is not of full rank.
404  // Here only rank-1 updates are implemented, to reduce the
405  // requirement for intermediate storage and improve accuracy
406  template<typename MatrixType, typename WDerived>
408  {
409  using numext::isfinite;
410  typedef typename MatrixType::Scalar Scalar;
411  typedef typename MatrixType::RealScalar RealScalar;
412 
413  const Index size = mat.rows();
414  eigen_assert(mat.cols() == size && w.size()==size);
415 
416  RealScalar alpha = 1;
417 
418  // Apply the update
419  for (Index j = 0; j < size; j++)
420  {
421  // Check for termination due to an original decomposition of low-rank
422  if (!(isfinite)(alpha))
423  break;
424 
425  // Update the diagonal terms
426  RealScalar dj = numext::real(mat.coeff(j,j));
427  Scalar wj = w.coeff(j);
428  RealScalar swj2 = sigma*numext::abs2(wj);
429  RealScalar gamma = dj*alpha + swj2;
430 
431  mat.coeffRef(j,j) += swj2/alpha;
432  alpha += swj2/dj;
433 
434 
435  // Update the terms of L
436  Index rs = size-j-1;
437  w.tail(rs) -= wj * mat.col(j).tail(rs);
438  if(gamma != 0)
439  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
440  }
441  return true;
442  }
443 
444  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
445  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
446  {
447  // Apply the permutation to the input w
448  tmp = transpositions * w;
449 
451  }
452 };
453 
454 template<> struct ldlt_inplace<Upper>
455 {
456  template<typename MatrixType, typename TranspositionType, typename Workspace>
457  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
458  {
459  Transpose<MatrixType> matt(mat);
460  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
461  }
462 
463  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
464  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
465  {
466  Transpose<MatrixType> matt(mat);
467  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
468  }
469 };
470 
471 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
472 {
475  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
476  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
477 };
478 
479 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
480 {
483  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
484  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
485 };
486 
487 } // end namespace internal
488 
491 template<typename MatrixType, int _UpLo>
492 template<typename InputType>
494 {
495  check_template_parameters();
496 
497  eigen_assert(a.rows()==a.cols());
498  const Index size = a.rows();
499 
500  m_matrix = a.derived();
501 
502  // Compute matrix L1 norm = max abs column sum.
503  m_l1_norm = RealScalar(0);
504  // TODO move this code to SelfAdjointView
505  for (Index col = 0; col < size; ++col) {
506  RealScalar abs_col_sum;
507  if (_UpLo == Lower)
508  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
509  else
510  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
511  if (abs_col_sum > m_l1_norm)
512  m_l1_norm = abs_col_sum;
513  }
514 
515  m_transpositions.resize(size);
516  m_isInitialized = false;
517  m_temporary.resize(size);
518  m_sign = internal::ZeroSign;
519 
520  m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
521 
522  m_isInitialized = true;
523  return *this;
524 }
525 
531 template<typename MatrixType, int _UpLo>
532 template<typename Derived>
534 {
535  typedef typename TranspositionType::StorageIndex IndexType;
536  const Index size = w.rows();
537  if (m_isInitialized)
538  {
539  eigen_assert(m_matrix.rows()==size);
540  }
541  else
542  {
543  m_matrix.resize(size,size);
544  m_matrix.setZero();
545  m_transpositions.resize(size);
546  for (Index i = 0; i < size; i++)
547  m_transpositions.coeffRef(i) = IndexType(i);
548  m_temporary.resize(size);
550  m_isInitialized = true;
551  }
552 
553  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
554 
555  return *this;
556 }
557 
558 #ifndef EIGEN_PARSED_BY_DOXYGEN
559 template<typename _MatrixType, int _UpLo>
560 template<typename RhsType, typename DstType>
561 void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
562 {
563  eigen_assert(rhs.rows() == rows());
564  // dst = P b
565  dst = m_transpositions * rhs;
566 
567  // dst = L^-1 (P b)
568  matrixL().solveInPlace(dst);
569 
570  // dst = D^-1 (L^-1 P b)
571  // more precisely, use pseudo-inverse of D (see bug 241)
572  using std::abs;
573  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
574  // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
575  // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
576  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
577  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
578  // diagonal element is not well justified and leads to numerical issues in some cases.
579  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
580  // Using numeric_limits::min() gives us more robustness to denormals.
582 
583  for (Index i = 0; i < vecD.size(); ++i)
584  {
585  if(abs(vecD(i)) > tolerance)
586  dst.row(i) /= vecD(i);
587  else
588  dst.row(i).setZero();
589  }
590 
591  // dst = L^-T (D^-1 L^-1 P b)
592  matrixU().solveInPlace(dst);
593 
594  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
595  dst = m_transpositions.transpose() * dst;
596 }
597 #endif
598 
612 template<typename MatrixType,int _UpLo>
613 template<typename Derived>
615 {
616  eigen_assert(m_isInitialized && "LDLT is not initialized.");
617  eigen_assert(m_matrix.rows() == bAndX.rows());
618 
619  bAndX = this->solve(bAndX);
620 
621  return true;
622 }
623 
627 template<typename MatrixType, int _UpLo>
629 {
630  eigen_assert(m_isInitialized && "LDLT is not initialized.");
631  const Index size = m_matrix.rows();
632  MatrixType res(size,size);
633 
634  // P
635  res.setIdentity();
636  res = transpositionsP() * res;
637  // L^* P
638  res = matrixU() * res;
639  // D(L^*P)
640  res = vectorD().real().asDiagonal() * res;
641  // L(DL^*P)
642  res = matrixL() * res;
643  // P^T (LDL^*P)
644  res = transpositionsP().transpose() * res;
645 
646  return res;
647 }
648 
653 template<typename MatrixType, unsigned int UpLo>
656 {
657  return LDLT<PlainObject,UpLo>(m_matrix);
658 }
659 
664 template<typename Derived>
667 {
668  return LDLT<PlainObject>(derived());
669 }
670 
671 } // end namespace Eigen
672 
673 #endif // EIGEN_LDLT_H
Matrix3f m
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:50
SCALAR Scalar
Definition: bench_gemm.cpp:33
static void check_template_parameters()
Definition: LDLT.h:267
#define EIGEN_STRONG_INLINE
Definition: Macros.h:494
TranspositionType m_transpositions
Definition: LDLT.h:280
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:666
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: LDLT.h:62
static bool update(MatrixType &mat, const TranspositionType &transpositions, Workspace &tmp, const WType &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:445
float real
Definition: datatypes.h:10
Scalar * b
Definition: benchVecAdd.cpp:17
LDLT & compute(const EigenBase< InputType > &matrix)
EIGEN_DEVICE_FUNC bool() isfinite(const T &x)
def update(text)
Definition: relicense.py:46
_MatrixType MatrixType
Definition: LDLT.h:53
#define min(a, b)
Definition: datatypes.h:19
Expression of the transpose of a matrix.
Definition: Transpose.h:52
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
Definition: LDLT.h:65
static const double sigma
Index rows() const
Definition: LDLT.h:245
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
LDLT & rankUpdate(const MatrixBase< Derived > &w, const RealScalar &alpha=1)
MatrixXf MatrixType
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
#define isfinite(X)
Definition: main.h:74
static bool updateInPlace(MatrixType &mat, MatrixBase< WDerived > &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:407
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
void setZero()
Definition: LDLT.h:135
RealScalar m_l1_norm
Definition: LDLT.h:279
MatrixType m_matrix
Definition: LDLT.h:278
Array33i a
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
internal::SignMatrix m_sign
Definition: LDLT.h:282
static EIGEN_STRONG_INLINE bool update(MatrixType &mat, TranspositionType &transpositions, Workspace &tmp, WType &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:464
MatrixType::StorageIndex StorageIndex
Definition: LDLT.h:64
LDLT()
Default Constructor.
Definition: LDLT.h:77
EIGEN_DEVICE_FUNC const SignReturnType sign() const
const mpreal gamma(const mpreal &x, mp_rnd_t r=mpreal::get_default_rnd())
Definition: mpreal.h:2262
bool solveInPlace(MatrixBase< Derived > &bAndX) const
Definition: LDLT.h:614
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
bool isNegative(void) const
Definition: LDLT.h:177
const TriangularView< const typename MatrixType::AdjointReturnType, UnitLower > MatrixL
Definition: LDLT.h:481
#define eigen_assert(x)
Definition: Macros.h:579
RealScalar alpha
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:184
RealScalar s
const TriangularView< const MatrixType, UnitLower > MatrixL
Definition: LDLT.h:473
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
Definition: LDLT.h:67
RealScalar rcond() const
Definition: LDLT.h:217
MatrixType reconstructedMatrix() const
Definition: LDLT.h:628
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:34
EIGEN_DEVICE_FUNC Index cols() const
Definition: EigenBase.h:62
static MatrixL getL(const MatrixType &m)
Definition: LDLT.h:475
RowVector3d w
const LDLT & adjoint() const
Definition: LDLT.h:243
internal::LDLT_Traits< MatrixType, UpLo > Traits
Definition: LDLT.h:70
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
static MatrixU getU(const MatrixType &m)
Definition: LDLT.h:476
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
static EIGEN_STRONG_INLINE bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
Definition: LDLT.h:457
Index cols() const
Definition: LDLT.h:246
TmpMatrixType m_temporary
Definition: LDLT.h:281
Eigen::Index Index
Definition: LDLT.h:63
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
ComputationInfo m_info
Definition: LDLT.h:284
EIGEN_DEVICE_FUNC SegmentReturnType tail(Index n)
Definition: DenseBase.h:950
bool isPositive() const
Definition: LDLT.h:170
static MatrixU getU(const MatrixType &m)
Definition: LDLT.h:484
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
DenseIndex ret
Definition: level1_impl.h:59
MatrixType::Scalar Scalar
Definition: LDLT.h:61
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
Expression of a triangular part in a matrix.
EIGEN_DEVICE_FUNC Index rows() const
Definition: EigenBase.h:59
m col(1)
const TriangularView< const MatrixType, UnitUpper > MatrixU
Definition: LDLT.h:482
const TriangularView< const typename MatrixType::AdjointReturnType, UnitUpper > MatrixU
Definition: LDLT.h:474
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
Pseudo expression representing a solving operation.
Definition: Solve.h:62
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
static bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
Definition: LDLT.h:294
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
#define abs(x)
Definition: datatypes.h:17
ComputationInfo
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:45
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationType
Definition: LDLT.h:68
bool m_isInitialized
Definition: LDLT.h:283
void swap(mpfr::mpreal &x, mpfr::mpreal &y)
Definition: mpreal.h:2986
std::ptrdiff_t j
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:655
static MatrixL getL(const MatrixType &m)
Definition: LDLT.h:483
ScalarWithExceptions conj(const ScalarWithExceptions &x)
Definition: exceptions.cpp:74


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:29