17#include "./InternalHeaderCheck.h"
22template <
typename MatrixType_,
int UpLo_>
23struct traits<LDLT<MatrixType_, UpLo_> > : traits<MatrixType_> {
24 typedef MatrixXpr XprKind;
25 typedef SolverStorage StorageKind;
26 typedef int StorageIndex;
30template <
typename MatrixType,
int UpLo>
34enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
62template <
typename MatrixType_,
int UpLo_>
65 typedef MatrixType_ MatrixType;
69 EIGEN_GENERIC_PUBLIC_INTERFACE(
LDLT)
71 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
80 typedef internal::LDLT_Traits<MatrixType, UpLo> Traits;
87 LDLT() : m_matrix(), m_transpositions(), m_sign(internal::ZeroSign), m_isInitialized(false) {}
97 m_transpositions(
size),
99 m_sign(internal::ZeroSign),
100 m_isInitialized(false) {}
108 template <
typename InputType>
110 : m_matrix(matrix.rows(), matrix.cols()),
111 m_transpositions(matrix.rows()),
112 m_temporary(matrix.rows()),
113 m_sign(internal::ZeroSign),
114 m_isInitialized(false) {
125 template <
typename InputType>
128 m_transpositions(matrix.rows()),
129 m_temporary(matrix.rows()),
130 m_sign(internal::ZeroSign),
131 m_isInitialized(false) {
141 inline typename Traits::MatrixU
matrixU()
const {
142 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
143 return Traits::getU(m_matrix);
147 inline typename Traits::MatrixL
matrixL()
const {
148 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
149 return Traits::getL(m_matrix);
155 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
156 return m_transpositions;
161 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
162 return m_matrix.diagonal();
167 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
168 return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
173 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
174 return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
177#ifdef EIGEN_PARSED_BY_DOXYGEN
193 template <
typename Rhs>
197 template <
typename Derived>
200 template <
typename InputType>
207 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
208 return internal::rcond_estimate_helper(m_l1_norm, *
this);
211 template <
typename Derived>
219 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
233 EIGEN_DEVICE_FUNC
inline EIGEN_CONSTEXPR
Index rows() const EIGEN_NOEXCEPT {
return m_matrix.rows(); }
234 EIGEN_DEVICE_FUNC
inline EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT {
return m_matrix.cols(); }
242 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
246#ifndef EIGEN_PARSED_BY_DOXYGEN
247 template <
typename RhsType,
typename DstType>
248 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
250 template <
bool Conjugate,
typename RhsType,
typename DstType>
251 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
255 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
264 RealScalar m_l1_norm;
265 TranspositionType m_transpositions;
266 TmpMatrixType m_temporary;
267 internal::SignMatrix m_sign;
268 bool m_isInitialized;
278struct ldlt_inplace<
Lower> {
279 template <
typename MatrixType,
typename TranspositionType,
typename Workspace>
280 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix&
sign) {
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;
291 transpositions.setIdentity();
294 else if (numext::real(mat.coeff(0, 0)) >
static_cast<RealScalar
>(0))
295 sign = PositiveSemiDef;
296 else if (numext::real(mat.coeff(0, 0)) <
static_cast<RealScalar
>(0))
297 sign = NegativeSemiDef;
303 for (Index k = 0; k < size; ++k) {
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;
309 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
310 if (k != index_of_biggest_in_corner) {
313 Index s = size - index_of_biggest_in_corner - 1;
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);
322 if (NumTraits<Scalar>::IsComplex)
323 mat.coeffRef(index_of_biggest_in_corner, k) = numext::conj(mat.coeff(index_of_biggest_in_corner, k));
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);
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);
345 RealScalar realAkk = numext::real(mat.coeffRef(k, k));
346 bool pivot_is_valid = (
abs(realAkk) > RealScalar(0));
348 if (k == 0 && !pivot_is_valid) {
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();
359 if ((rs > 0) && pivot_is_valid)
362 ret = ret && (A21.array() == Scalar(0)).all();
364 if (found_zero_pivot && pivot_is_valid)
366 else if (!pivot_is_valid)
367 found_zero_pivot =
true;
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))
375 sign = PositiveSemiDef;
376 else if (realAkk <
static_cast<RealScalar
>(0))
377 sign = NegativeSemiDef;
391 template <
typename MatrixType,
typename WDerived>
392 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w,
393 const typename MatrixType::RealScalar& sigma = 1) {
394 using numext::isfinite;
395 typedef typename MatrixType::Scalar Scalar;
396 typedef typename MatrixType::RealScalar RealScalar;
398 const Index size = mat.rows();
399 eigen_assert(mat.cols() == size && w.size() == size);
401 RealScalar alpha = 1;
404 for (Index j = 0; j < size; j++) {
409 RealScalar dj = numext::real(mat.coeff(j, j));
410 Scalar wj = w.coeff(j);
411 RealScalar swj2 = sigma * numext::abs2(wj);
412 RealScalar gamma = dj * alpha + swj2;
414 mat.coeffRef(j, j) += swj2 / alpha;
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);
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) {
429 tmp = transpositions * w;
431 return ldlt_inplace<Lower>::updateInPlace(mat, tmp, sigma);
436struct ldlt_inplace<
Upper> {
437 template <
typename MatrixType,
typename TranspositionType,
typename Workspace>
438 static EIGEN_STRONG_INLINE
bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp,
440 Transpose<MatrixType> matt(mat);
441 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp,
sign);
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) {
447 Transpose<MatrixType> matt(mat);
448 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
452template <
typename MatrixType>
453struct LDLT_Traits<MatrixType,
Lower> {
454 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
455 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
456 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m); }
457 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m.adjoint()); }
460template <
typename MatrixType>
461struct LDLT_Traits<MatrixType,
Upper> {
462 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
463 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
464 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m.adjoint()); }
465 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m); }
472template <
typename MatrixType,
int UpLo_>
473template <
typename InputType>
481 m_l1_norm = RealScalar(0);
483 for (
Index col = 0; col < size; ++col) {
484 RealScalar abs_col_sum;
487 m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
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;
494 m_transpositions.resize(size);
495 m_isInitialized =
false;
496 m_temporary.resize(size);
497 m_sign = internal::ZeroSign;
499 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ?
Success
502 m_isInitialized =
true;
511template <
typename MatrixType,
int UpLo_>
512template <
typename Derived>
515 typedef typename TranspositionType::StorageIndex IndexType;
517 if (m_isInitialized) {
518 eigen_assert(m_matrix.rows() == size);
520 m_matrix.resize(size, size);
522 m_transpositions.resize(size);
523 for (
Index i = 0; i < size; i++) m_transpositions.coeffRef(i) = IndexType(i);
524 m_temporary.resize(size);
525 m_sign = sigma >= 0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
526 m_isInitialized =
true;
529 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
534#ifndef EIGEN_PARSED_BY_DOXYGEN
535template <
typename MatrixType_,
int UpLo_>
536template <
typename RhsType,
typename DstType>
538 _solve_impl_transposed<true>(rhs, dst);
541template <
typename MatrixType_,
int UpLo_>
542template <
bool Conjugate,
typename RhsType,
typename DstType>
543void LDLT<MatrixType_, UpLo_>::_solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const {
545 dst = m_transpositions * rhs;
549 matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
555 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
562 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
563 for (Index i = 0; i < vecD.size(); ++i) {
564 if (abs(vecD(i)) > tolerance)
565 dst.row(i) /= vecD(i);
567 dst.row(i).setZero();
572 matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
576 dst = m_transpositions.transpose() * dst;
593template <
typename MatrixType,
int UpLo_>
594template <
typename Derived>
595bool LDLT<MatrixType, UpLo_>::solveInPlace(MatrixBase<Derived>& bAndX)
const {
596 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
597 eigen_assert(m_matrix.rows() == bAndX.rows());
599 bAndX = this->solve(bAndX);
607template <
typename MatrixType,
int UpLo_>
609 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
610 const Index size = m_matrix.rows();
611 MatrixType res(size, size);
615 res = transpositionsP() * res;
617 res = matrixU() * res;
619 res = vectorD().real().asDiagonal() * res;
621 res = matrixL() * res;
623 res = transpositionsP().transpose() * res;
632template <
typename MatrixType,
unsigned int UpLo>
642template <
typename Derived>
internal::traits< Homogeneous< MatrixType, Direction_ > >::Scalar Scalar
Definition DenseBase.h:62
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition EigenBase.h:59
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
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
void setZero()
Definition LDLT.h:138
const MatrixType & matrixLDLT() const
Definition LDLT.h:218
RealScalar rcond() const
Definition LDLT.h:206
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
bool isNegative(void) const
Definition LDLT.h:172
LDLT()
Default Constructor.
Definition LDLT.h:87
bool isPositive() const
Definition LDLT.h:166
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition LDLT.h:109
LDLT(Index size)
Default Constructor with memory preallocation.
Definition LDLT.h:95
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition LDLT.h:126
MatrixType reconstructedMatrix() const
Definition LDLT.h:608
ComputationInfo info() const
Reports whether previous computation was successful.
Definition LDLT.h:241
const LDLT & adjoint() const
Definition LDLT.h:231
Traits::MatrixU matrixU() const
Definition LDLT.h:141
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Permutation matrix.
Definition PermutationMatrix.h:280
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition SelfAdjointView.h:51
Pseudo expression representing a solving operation.
Definition Solve.h:62
A base class for matrix decomposition and solvers.
Definition SolverBase.h:72
LDLT< MatrixType_, UpLo_ > & derived()
Definition EigenBase.h:49
Represents a sequence of transpositions (row/column interchange)
Definition Transpositions.h:141
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
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition EigenBase.h:61
Derived & derived()
Definition EigenBase.h:49
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition EigenBase.h:59
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition EigenBase.h:64