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 // IWYU pragma: private
17 #include "./InternalHeaderCheck.h"
18 
19 namespace Eigen {
20 
21 namespace internal {
22 template <typename MatrixType_, int UpLo_>
23 struct traits<LDLT<MatrixType_, UpLo_> > : traits<MatrixType_> {
24  typedef MatrixXpr XprKind;
26  typedef int StorageIndex;
27  enum { Flags = 0 };
28 };
29 
30 template <typename MatrixType, int UpLo>
31 struct LDLT_Traits;
32 
33 // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
35 } // namespace internal
36 
62 template <typename MatrixType_, int UpLo_>
63 class LDLT : public SolverBase<LDLT<MatrixType_, UpLo_> > {
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 
81 
88 
95  explicit LDLT(Index size)
96  : m_matrix(size, size),
100  m_isInitialized(false) {}
101 
108  template <typename InputType>
110  : m_matrix(matrix.rows(), matrix.cols()),
114  m_isInitialized(false) {
115  compute(matrix.derived());
116  }
117 
125  template <typename InputType>
127  : m_matrix(matrix.derived()),
131  m_isInitialized(false) {
132  compute(matrix.derived());
133  }
134 
138  void setZero() { m_isInitialized = false; }
139 
141  inline typename Traits::MatrixU matrixU() const {
142  eigen_assert(m_isInitialized && "LDLT is not initialized.");
143  return Traits::getU(m_matrix);
144  }
145 
147  inline typename Traits::MatrixL matrixL() const {
148  eigen_assert(m_isInitialized && "LDLT is not initialized.");
149  return Traits::getL(m_matrix);
150  }
151 
154  inline const TranspositionType& transpositionsP() const {
155  eigen_assert(m_isInitialized && "LDLT is not initialized.");
156  return m_transpositions;
157  }
158 
161  eigen_assert(m_isInitialized && "LDLT is not initialized.");
162  return m_matrix.diagonal();
163  }
164 
166  inline bool isPositive() const {
167  eigen_assert(m_isInitialized && "LDLT is not initialized.");
169  }
170 
172  inline bool isNegative(void) const {
173  eigen_assert(m_isInitialized && "LDLT is not initialized.");
175  }
176 
177 #ifdef EIGEN_PARSED_BY_DOXYGEN
193  template <typename Rhs>
194  inline const Solve<LDLT, Rhs> solve(const MatrixBase<Rhs>& b) const;
195 #endif
196 
197  template <typename Derived>
198  bool solveInPlace(MatrixBase<Derived>& bAndX) const;
199 
200  template <typename InputType>
202 
206  RealScalar rcond() const {
207  eigen_assert(m_isInitialized && "LDLT is not initialized.");
209  }
210 
211  template <typename Derived>
213 
218  inline const MatrixType& matrixLDLT() const {
219  eigen_assert(m_isInitialized && "LDLT is not initialized.");
220  return m_matrix;
221  }
222 
224 
231  const LDLT& adjoint() const { return *this; }
232 
235 
242  eigen_assert(m_isInitialized && "LDLT is not initialized.");
243  return m_info;
244  }
245 
246 #ifndef EIGEN_PARSED_BY_DOXYGEN
247  template <typename RhsType, typename DstType>
248  void _solve_impl(const RhsType& rhs, DstType& dst) const;
249 
250  template <bool Conjugate, typename RhsType, typename DstType>
251  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
252 #endif
253 
254  protected:
256 
257 
270 };
271 
272 namespace internal {
273 
274 template <int UpLo>
276 
277 template <>
279  template <typename MatrixType, typename TranspositionType, typename Workspace>
280  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign) {
281  using std::abs;
282  typedef typename MatrixType::Scalar Scalar;
283  typedef typename MatrixType::RealScalar RealScalar;
284  typedef typename TranspositionType::StorageIndex IndexType;
285  eigen_assert(mat.rows() == mat.cols());
286  const Index size = mat.rows();
287  bool found_zero_pivot = false;
288  bool ret = true;
289 
290  if (size <= 1) {
291  transpositions.setIdentity();
292  if (size == 0)
293  sign = ZeroSign;
294  else if (numext::real(mat.coeff(0, 0)) > static_cast<RealScalar>(0))
296  else if (numext::real(mat.coeff(0, 0)) < static_cast<RealScalar>(0))
298  else
299  sign = ZeroSign;
300  return true;
301  }
302 
303  for (Index k = 0; k < size; ++k) {
304  // Find largest diagonal element
305  Index index_of_biggest_in_corner;
306  mat.diagonal().tail(size - k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
307  index_of_biggest_in_corner += k;
308 
309  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
310  if (k != index_of_biggest_in_corner) {
311  // apply the transposition while taking care to consider only
312  // the lower triangular part
313  Index s = size - index_of_biggest_in_corner - 1; // trailing size after the biggest element
314  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
315  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
316  std::swap(mat.coeffRef(k, k), mat.coeffRef(index_of_biggest_in_corner, index_of_biggest_in_corner));
317  for (Index i = k + 1; i < index_of_biggest_in_corner; ++i) {
318  Scalar tmp = mat.coeffRef(i, k);
319  mat.coeffRef(i, k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner, i));
320  mat.coeffRef(index_of_biggest_in_corner, i) = numext::conj(tmp);
321  }
323  mat.coeffRef(index_of_biggest_in_corner, k) = numext::conj(mat.coeff(index_of_biggest_in_corner, k));
324  }
325 
326  // partition the matrix:
327  // A00 | - | -
328  // lu = A10 | A11 | -
329  // A20 | A21 | A22
330  Index rs = size - k - 1;
331  Block<MatrixType, Dynamic, 1> A21(mat, k + 1, k, rs, 1);
332  Block<MatrixType, 1, Dynamic> A10(mat, k, 0, 1, k);
333  Block<MatrixType, Dynamic, Dynamic> A20(mat, k + 1, 0, rs, k);
334 
335  if (k > 0) {
336  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
337  mat.coeffRef(k, k) -= (A10 * temp.head(k)).value();
338  if (rs > 0) A21.noalias() -= A20 * temp.head(k);
339  }
340 
341  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
342  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
343  // we should only make sure that we do not introduce INF or NaN values.
344  // Remark that LAPACK also uses 0 as the cutoff value.
345  RealScalar realAkk = numext::real(mat.coeffRef(k, k));
346  bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
347 
348  if (k == 0 && !pivot_is_valid) {
349  // The entire diagonal is zero, there is nothing more to do
350  // except filling the transpositions, and checking whether the matrix is zero.
351  sign = ZeroSign;
352  for (Index j = 0; j < size; ++j) {
353  transpositions.coeffRef(j) = IndexType(j);
354  ret = ret && (mat.col(j).tail(size - j - 1).array() == Scalar(0)).all();
355  }
356  return ret;
357  }
358 
359  if ((rs > 0) && pivot_is_valid)
360  A21 /= realAkk;
361  else if (rs > 0)
362  ret = ret && (A21.array() == Scalar(0)).all();
363 
364  if (found_zero_pivot && pivot_is_valid)
365  ret = false; // factorization failed
366  else if (!pivot_is_valid)
367  found_zero_pivot = true;
368 
369  if (sign == PositiveSemiDef) {
370  if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
371  } else if (sign == NegativeSemiDef) {
372  if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
373  } else if (sign == ZeroSign) {
374  if (realAkk > static_cast<RealScalar>(0))
376  else if (realAkk < static_cast<RealScalar>(0))
378  }
379  }
380 
381  return ret;
382  }
383 
384  // Reference for the algorithm: Davis and Hager, "Multiple Rank
385  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
386  // Trivial rearrangements of their computations (Timothy E. Holy)
387  // allow their algorithm to work for rank-1 updates even if the
388  // original matrix is not of full rank.
389  // Here only rank-1 updates are implemented, to reduce the
390  // requirement for intermediate storage and improve accuracy
391  template <typename MatrixType, typename WDerived>
393  const typename MatrixType::RealScalar& sigma = 1) {
394  using numext::isfinite;
395  typedef typename MatrixType::Scalar Scalar;
396  typedef typename MatrixType::RealScalar RealScalar;
397 
398  const Index size = mat.rows();
399  eigen_assert(mat.cols() == size && w.size() == size);
400 
401  RealScalar alpha = 1;
402 
403  // Apply the update
404  for (Index j = 0; j < size; j++) {
405  // Check for termination due to an original decomposition of low-rank
406  if (!(isfinite)(alpha)) break;
407 
408  // Update the diagonal terms
410  Scalar wj = w.coeff(j);
411  RealScalar swj2 = sigma * numext::abs2(wj);
412  RealScalar gamma = dj * alpha + swj2;
413 
414  mat.coeffRef(j, j) += swj2 / alpha;
415  alpha += swj2 / dj;
416 
417  // Update the terms of L
418  Index rs = size - j - 1;
419  w.tail(rs) -= wj * mat.col(j).tail(rs);
420  if (!numext::is_exactly_zero(gamma)) mat.col(j).tail(rs) += (sigma * numext::conj(wj) / gamma) * w.tail(rs);
421  }
422  return true;
423  }
424 
425  template <typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
426  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w,
427  const typename MatrixType::RealScalar& sigma = 1) {
428  // Apply the permutation to the input w
429  tmp = transpositions * w;
430 
432  }
433 };
434 
435 template <>
437  template <typename MatrixType, typename TranspositionType, typename Workspace>
438  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp,
439  SignMatrix& sign) {
441  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
442  }
443 
444  template <typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
445  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w,
446  const typename MatrixType::RealScalar& sigma = 1) {
448  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
449  }
450 };
451 
452 template <typename MatrixType>
456  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
457  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
458 };
459 
460 template <typename MatrixType>
464  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
465  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
466 };
467 
468 } // end namespace internal
469 
472 template <typename MatrixType, int UpLo_>
473 template <typename InputType>
475  eigen_assert(a.rows() == a.cols());
476  const Index size = a.rows();
477 
478  m_matrix = a.derived();
479 
480  // Compute matrix L1 norm = max abs column sum.
481  m_l1_norm = RealScalar(0);
482  // TODO move this code to SelfAdjointView
483  for (Index col = 0; col < size; ++col) {
484  RealScalar abs_col_sum;
485  if (UpLo_ == Lower)
486  abs_col_sum =
487  m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
488  else
489  abs_col_sum =
490  m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
491  if (abs_col_sum > m_l1_norm) m_l1_norm = abs_col_sum;
492  }
493 
495  m_isInitialized = false;
498 
500  : NumericalIssue;
501 
502  m_isInitialized = true;
503  return *this;
504 }
505 
511 template <typename MatrixType, int UpLo_>
512 template <typename Derived>
515  typedef typename TranspositionType::StorageIndex IndexType;
516  const Index size = w.rows();
517  if (m_isInitialized) {
518  eigen_assert(m_matrix.rows() == size);
519  } else {
520  m_matrix.resize(size, size);
521  m_matrix.setZero();
523  for (Index i = 0; i < size; i++) m_transpositions.coeffRef(i) = IndexType(i);
526  m_isInitialized = true;
527  }
528 
530 
531  return *this;
532 }
533 
534 #ifndef EIGEN_PARSED_BY_DOXYGEN
535 template <typename MatrixType_, int UpLo_>
536 template <typename RhsType, typename DstType>
537 void LDLT<MatrixType_, UpLo_>::_solve_impl(const RhsType& rhs, DstType& dst) const {
538  _solve_impl_transposed<true>(rhs, dst);
539 }
540 
541 template <typename MatrixType_, int UpLo_>
542 template <bool Conjugate, typename RhsType, typename DstType>
543 void LDLT<MatrixType_, UpLo_>::_solve_impl_transposed(const RhsType& rhs, DstType& dst) const {
544  // dst = P b
545  dst = m_transpositions * rhs;
546 
547  // dst = L^-1 (P b)
548  // dst = L^-*T (P b)
549  matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
550 
551  // dst = D^-* (L^-1 P b)
552  // dst = D^-1 (L^-*T P b)
553  // more precisely, use pseudo-inverse of D (see bug 241)
554  using std::abs;
556  // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
557  // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
558  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1)
559  // / NumTraits<RealScalar>::highest()); However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the
560  // highest diagonal element is not well justified and leads to numerical issues in some cases. Moreover, Lapack's
561  // xSYTRS routines use 0 for the tolerance. Using numeric_limits::min() gives us more robustness to denormals.
563  for (Index i = 0; i < vecD.size(); ++i) {
564  if (abs(vecD(i)) > tolerance)
565  dst.row(i) /= vecD(i);
566  else
567  dst.row(i).setZero();
568  }
569 
570  // dst = L^-* (D^-* L^-1 P b)
571  // dst = L^-T (D^-1 L^-*T P b)
572  matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
573 
574  // dst = P^T (L^-* D^-* L^-1 P b) = A^-1 b
575  // dst = P^-T (L^-T D^-1 L^-*T P b) = A^-1 b
576  dst = m_transpositions.transpose() * dst;
577 }
578 #endif
579 
593 template <typename MatrixType, int UpLo_>
594 template <typename Derived>
596  eigen_assert(m_isInitialized && "LDLT is not initialized.");
597  eigen_assert(m_matrix.rows() == bAndX.rows());
598 
599  bAndX = this->solve(bAndX);
600 
601  return true;
602 }
603 
607 template <typename MatrixType, int UpLo_>
609  eigen_assert(m_isInitialized && "LDLT is not initialized.");
610  const Index size = m_matrix.rows();
612 
613  // P
614  res.setIdentity();
615  res = transpositionsP() * res;
616  // L^* P
617  res = matrixU() * res;
618  // D(L^*P)
619  res = vectorD().real().asDiagonal() * res;
620  // L(DL^*P)
621  res = matrixL() * res;
622  // P^T (LDL^*P)
624 
625  return res;
626 }
627 
632 template <typename MatrixType, unsigned int UpLo>
636 }
637 
642 template <typename Derived>
644  return LDLT<PlainObject>(derived());
645 }
646 
647 } // end namespace Eigen
648 
649 #endif // EIGEN_LDLT_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
Definition: Macros.h:1149
#define EIGEN_NOEXCEPT
Definition: Macros.h:1267
#define EIGEN_CONSTEXPR
Definition: Macros.h:758
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:892
#define eigen_assert(x)
Definition: Macros.h:910
#define EIGEN_STRONG_INLINE
Definition: Macros.h:834
m col(1)
RowVector3d w
Definition: Matrix_resize_int.cpp:3
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:74
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:110
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:68
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:154
Traits::MatrixL matrixL() const
Definition: LDLT.h:147
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:160
ComputationInfo m_info
Definition: LDLT.h:269
internal::SignMatrix m_sign
Definition: LDLT.h:267
const MatrixType & matrixLDLT() const
Definition: LDLT.h:218
bool m_isInitialized
Definition: LDLT.h:268
void setZero()
Definition: LDLT.h:138
void _solve_impl(const RhsType &rhs, DstType &dst) const
Definition: LDLT.h:537
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: LDLT.h:234
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: LDLT.h:233
TranspositionType m_transpositions
Definition: LDLT.h:265
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
Definition: LDLT.h:77
RealScalar rcond() const
Definition: LDLT.h:206
internal::LDLT_Traits< MatrixType, UpLo > Traits
Definition: LDLT.h:80
MatrixType m_matrix
Definition: LDLT.h:263
bool isNegative(void) const
Definition: LDLT.h:172
LDLT()
Default Constructor.
Definition: LDLT.h:87
bool isPositive() const
Definition: LDLT.h:166
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
Definition: LDLT.h:543
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:109
SolverBase< LDLT > Base
Definition: LDLT.h:66
RealScalar m_l1_norm
Definition: LDLT.h:264
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:95
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:126
Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
Definition: LDLT.h:75
TmpMatrixType m_temporary
Definition: LDLT.h:266
MatrixType reconstructedMatrix() const
Definition: LDLT.h:608
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:241
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
bool solveInPlace(MatrixBase< Derived > &bAndX) const
Definition: LDLT.h:595
const LDLT & adjoint() const
Definition: LDLT.h:231
@ MaxColsAtCompileTime
Definition: LDLT.h:72
@ MaxRowsAtCompileTime
Definition: LDLT.h:71
@ UpLo
Definition: LDLT.h:73
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:643
Permutation matrix.
Definition: PermutationMatrix.h:280
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:634
Pseudo expression representing a solving operation.
Definition: Solve.h:62
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:72
internal::traits< LDLT< MatrixType_, UpLo_ > >::Scalar Scalar
Definition: SolverBase.h:75
constexpr EIGEN_DEVICE_FUNC LDLT< MatrixType_, UpLo_ > & derived()
Definition: EigenBase.h:49
const Solve< LDLT< MatrixType_, UpLo_ >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SolverBase.h:106
Scalar coeff(Index row, Index col) const
Definition: SparseMatrix.h:211
const ConstDiagonalReturnType diagonal() const
Definition: SparseMatrix.h:757
void swap(SparseMatrix &other)
Definition: SparseMatrix.h:829
Index cols() const
Definition: SparseMatrix.h:161
Scalar & coeffRef(Index row, Index col)
Definition: SparseMatrix.h:275
Index rows() const
Definition: SparseMatrix.h:159
Expression of the transpose of a matrix.
Definition: Transpose.h:56
Transpose< TranspositionsBase > transpose() const
Definition: Transpositions.h:95
void resize(Index newSize)
Definition: Transpositions.h:63
StorageIndex & coeffRef(Index i)
Definition: Transpositions.h:47
IndicesType::Scalar StorageIndex
Definition: Transpositions.h:147
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:167
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
float real
Definition: datatypes.h:10
#define min(a, b)
Definition: datatypes.h:22
static constexpr Eigen::internal::all_t all
Definition: IndexedViewHelper.h:86
ComputationInfo
Definition: Constants.h:438
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
RealScalar s
Definition: level1_cplx_impl.h:130
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
RealScalar alpha
Definition: level1_cplx_impl.h:151
EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:117
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
#define isfinite(X)
Definition: main.h:111
SignMatrix
Definition: LDLT.h:34
@ PositiveSemiDef
Definition: LDLT.h:34
@ ZeroSign
Definition: LDLT.h:34
@ NegativeSemiDef
Definition: LDLT.h:34
@ Indefinite
Definition: LDLT.h:34
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
Definition: ConditionEstimator.h:157
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE bool() isfinite(const Eigen::bfloat16 &h)
Definition: BFloat16.h:752
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
Extend namespace for flags.
Definition: fsi_chan_precond_driver.cc:56
T sign(T x)
Definition: cxx11_tensor_builtins_sycl.cpp:172
int sigma
Definition: calibrate.py:179
Definition: Eigen_Colamd.h:49
Mdouble gamma(Mdouble gamma_in)
This is the gamma function returns the true value for the half integer value.
Definition: ExtendedMath.cc:116
Definition: EigenBase.h:33
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:64
Definition: Constants.h:534
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: Constants.h:525
static MatrixL getL(const MatrixType &m)
Definition: LDLT.h:456
const TriangularView< const typename MatrixType::AdjointReturnType, UnitUpper > MatrixU
Definition: LDLT.h:455
static MatrixU getU(const MatrixType &m)
Definition: LDLT.h:457
const TriangularView< const MatrixType, UnitLower > MatrixL
Definition: LDLT.h:454
const TriangularView< const typename MatrixType::AdjointReturnType, UnitLower > MatrixL
Definition: LDLT.h:462
static MatrixL getL(const MatrixType &m)
Definition: LDLT.h:464
const TriangularView< const MatrixType, UnitUpper > MatrixU
Definition: LDLT.h:463
static MatrixU getU(const MatrixType &m)
Definition: LDLT.h:465
Definition: LDLT.h:31
static bool updateInPlace(MatrixType &mat, MatrixBase< WDerived > &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:392
static bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
Definition: LDLT.h:280
static bool update(MatrixType &mat, const TranspositionType &transpositions, Workspace &tmp, const WType &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:426
static EIGEN_STRONG_INLINE bool update(MatrixType &mat, TranspositionType &transpositions, Workspace &tmp, WType &w, const typename MatrixType::RealScalar &sigma=1)
Definition: LDLT.h:445
static EIGEN_STRONG_INLINE bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
Definition: LDLT.h:438
Definition: LDLT.h:275
SolverStorage StorageKind
Definition: LDLT.h:25
Definition: ForwardDeclarations.h:21
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2