19 template<
typename _MatrixType,
int _UpLo>
struct traits<LDLT<_MatrixType, _UpLo> >
22 typedef MatrixXpr XprKind;
23 typedef SolverStorage StorageKind;
24 typedef int StorageIndex;
28 template<
typename MatrixType,
int UpLo>
struct LDLT_Traits;
31 enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
59template<
typename _MatrixType,
int _UpLo>
class LDLT
63 typedef _MatrixType MatrixType;
67 EIGEN_GENERIC_PUBLIC_INTERFACE(
LDLT)
69 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
70 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
78 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
88 m_sign(internal::ZeroSign),
89 m_isInitialized(false)
100 m_transpositions(
size),
102 m_sign(internal::ZeroSign),
103 m_isInitialized(false)
112 template<
typename InputType>
114 : m_matrix(matrix.rows(), matrix.cols()),
115 m_transpositions(matrix.rows()),
116 m_temporary(matrix.rows()),
117 m_sign(internal::ZeroSign),
118 m_isInitialized(false)
129 template<
typename InputType>
132 m_transpositions(matrix.rows()),
133 m_temporary(matrix.rows()),
134 m_sign(internal::ZeroSign),
135 m_isInitialized(false)
145 m_isInitialized =
false;
149 inline typename Traits::MatrixU
matrixU()
const
151 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
152 return Traits::getU(m_matrix);
156 inline typename Traits::MatrixL
matrixL()
const
158 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
159 return Traits::getL(m_matrix);
166 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
167 return m_transpositions;
173 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
174 return m_matrix.diagonal();
180 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
181 return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
187 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
188 return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
191 #ifdef EIGEN_PARSED_BY_DOXYGEN
207 template<
typename Rhs>
212 template<
typename Derived>
215 template<
typename InputType>
223 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
224 return internal::rcond_estimate_helper(m_l1_norm, *
this);
227 template <
typename Derived>
236 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
249 EIGEN_DEVICE_FUNC
inline EIGEN_CONSTEXPR
Index rows() const EIGEN_NOEXCEPT {
return m_matrix.rows(); }
250 EIGEN_DEVICE_FUNC
inline EIGEN_CONSTEXPR
Index cols() const EIGEN_NOEXCEPT {
return m_matrix.cols(); }
259 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
263 #ifndef EIGEN_PARSED_BY_DOXYGEN
264 template<
typename RhsType,
typename DstType>
265 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
267 template<
bool Conjugate,
typename RhsType,
typename DstType>
268 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
273 static void check_template_parameters()
275 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
285 RealScalar m_l1_norm;
286 TranspositionType m_transpositions;
287 TmpMatrixType m_temporary;
288 internal::SignMatrix m_sign;
289 bool m_isInitialized;
295template<
int UpLo>
struct ldlt_inplace;
297template<>
struct ldlt_inplace<
Lower>
299 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
300 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
303 typedef typename MatrixType::Scalar Scalar;
304 typedef typename MatrixType::RealScalar RealScalar;
305 typedef typename TranspositionType::StorageIndex IndexType;
306 eigen_assert(mat.rows()==mat.cols());
307 const Index size = mat.rows();
308 bool found_zero_pivot =
false;
313 transpositions.setIdentity();
314 if(size==0)
sign = ZeroSign;
315 else if (numext::real(mat.coeff(0,0)) >
static_cast<RealScalar
>(0) )
sign = PositiveSemiDef;
316 else if (numext::real(mat.coeff(0,0)) <
static_cast<RealScalar
>(0))
sign = NegativeSemiDef;
317 else sign = ZeroSign;
321 for (
Index k = 0; k < size; ++k)
324 Index index_of_biggest_in_corner;
325 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
326 index_of_biggest_in_corner += k;
328 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
329 if(k != index_of_biggest_in_corner)
333 Index s = size-index_of_biggest_in_corner-1;
334 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
335 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
336 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
337 for(
Index i=k+1;i<index_of_biggest_in_corner;++i)
339 Scalar tmp = mat.coeffRef(i,k);
340 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
341 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
343 if(NumTraits<Scalar>::IsComplex)
344 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
351 Index rs = size - k - 1;
352 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
353 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
354 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
358 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
359 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
361 A21.noalias() -= A20 * temp.head(k);
368 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
369 bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
371 if(k==0 && !pivot_is_valid)
376 for(
Index j = 0; j<size; ++j)
378 transpositions.coeffRef(j) = IndexType(j);
379 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).
all();
384 if((rs>0) && pivot_is_valid)
387 ret = ret && (A21.array()==Scalar(0)).
all();
389 if(found_zero_pivot && pivot_is_valid) ret =
false;
390 else if(!pivot_is_valid) found_zero_pivot =
true;
392 if (sign == PositiveSemiDef) {
393 if (realAkk <
static_cast<RealScalar
>(0))
sign = Indefinite;
394 }
else if (sign == NegativeSemiDef) {
395 if (realAkk >
static_cast<RealScalar
>(0))
sign = Indefinite;
396 }
else if (sign == ZeroSign) {
397 if (realAkk >
static_cast<RealScalar
>(0))
sign = PositiveSemiDef;
398 else if (realAkk <
static_cast<RealScalar
>(0))
sign = NegativeSemiDef;
412 template<
typename MatrixType,
typename WDerived>
413 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w,
const typename MatrixType::RealScalar& sigma=1)
415 using numext::isfinite;
416 typedef typename MatrixType::Scalar Scalar;
417 typedef typename MatrixType::RealScalar RealScalar;
419 const Index size = mat.rows();
420 eigen_assert(mat.cols() == size && w.size()==size);
422 RealScalar alpha = 1;
425 for (
Index j = 0; j < size; j++)
432 RealScalar dj = numext::real(mat.coeff(j,j));
433 Scalar wj = w.coeff(j);
434 RealScalar swj2 = sigma*numext::abs2(wj);
435 RealScalar gamma = dj*alpha + swj2;
437 mat.coeffRef(j,j) += swj2/alpha;
443 w.tail(rs) -= wj * mat.col(j).tail(rs);
445 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
450 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
451 static bool update(MatrixType& mat,
const TranspositionType& transpositions, Workspace& tmp,
const WType& w,
const typename MatrixType::RealScalar& sigma=1)
454 tmp = transpositions * w;
456 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
460template<>
struct ldlt_inplace<
Upper>
462 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
463 static EIGEN_STRONG_INLINE
bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
465 Transpose<MatrixType> matt(mat);
466 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
469 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
470 static EIGEN_STRONG_INLINE
bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w,
const typename MatrixType::RealScalar& sigma=1)
472 Transpose<MatrixType> matt(mat);
473 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
477template<
typename MatrixType>
struct LDLT_Traits<MatrixType,
Lower>
479 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
480 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
481 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m); }
482 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m.adjoint()); }
485template<
typename MatrixType>
struct LDLT_Traits<MatrixType,
Upper>
487 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
488 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
489 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m.adjoint()); }
490 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m); }
497template<
typename MatrixType,
int _UpLo>
498template<
typename InputType>
501 check_template_parameters();
509 m_l1_norm = RealScalar(0);
511 for (
Index col = 0; col < size; ++col) {
512 RealScalar abs_col_sum;
514 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
516 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
517 if (abs_col_sum > m_l1_norm)
518 m_l1_norm = abs_col_sum;
521 m_transpositions.resize(size);
522 m_isInitialized =
false;
523 m_temporary.resize(size);
524 m_sign = internal::ZeroSign;
526 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ?
Success :
NumericalIssue;
528 m_isInitialized =
true;
537template<
typename MatrixType,
int _UpLo>
538template<
typename Derived>
541 typedef typename TranspositionType::StorageIndex IndexType;
545 eigen_assert(m_matrix.rows()==size);
549 m_matrix.resize(size,size);
551 m_transpositions.resize(size);
552 for (
Index i = 0; i < size; i++)
553 m_transpositions.coeffRef(i) = IndexType(i);
554 m_temporary.resize(size);
555 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
556 m_isInitialized =
true;
559 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
564#ifndef EIGEN_PARSED_BY_DOXYGEN
565template<
typename _MatrixType,
int _UpLo>
566template<
typename RhsType,
typename DstType>
569 _solve_impl_transposed<true>(rhs, dst);
572template<
typename _MatrixType,
int _UpLo>
573template<
bool Conjugate,
typename RhsType,
typename DstType>
574void LDLT<_MatrixType,_UpLo>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const
577 dst = m_transpositions * rhs;
581 matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
587 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
595 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
596 for (Index i = 0; i < vecD.size(); ++i)
598 if(abs(vecD(i)) > tolerance)
599 dst.row(i) /= vecD(i);
601 dst.row(i).setZero();
606 matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
610 dst = m_transpositions.transpose() * dst;
627template<
typename MatrixType,
int _UpLo>
628template<
typename Derived>
629bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX)
const
631 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
632 eigen_assert(m_matrix.rows() == bAndX.rows());
634 bAndX = this->solve(bAndX);
642template<
typename MatrixType,
int _UpLo>
645 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
646 const Index size = m_matrix.rows();
647 MatrixType res(size,size);
651 res = transpositionsP() * res;
653 res = matrixU() * res;
655 res = vectorD().real().asDiagonal() * res;
657 res = matrixL() * res;
659 res = transpositionsP().transpose() * res;
668template<
typename MatrixType,
unsigned int UpLo>
679template<
typename Derived>
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:60
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:65
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:61
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:98
LDLT()
Default Constructor.
Definition: LDLT.h:85
Traits::MatrixU matrixU() const
Definition: LDLT.h:149
bool isPositive() const
Definition: LDLT.h:178
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:257
void setZero()
Definition: LDLT.h:143
const MatrixType & matrixLDLT() const
Definition: LDLT.h:234
bool isNegative(void) const
Definition: LDLT.h:185
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:171
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:113
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:130
const LDLT & adjoint() const
Definition: LDLT.h:247
MatrixType reconstructedMatrix() const
Definition: LDLT.h:643
RealScalar rcond() const
Definition: LDLT.h:221
Traits::MatrixL matrixL() const
Definition: LDLT.h:156
const TranspositionType & transpositionsP() const
Definition: LDLT.h:164
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Permutation matrix.
Definition: PermutationMatrix.h:298
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:63
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:69
LDLT< _MatrixType, _UpLo > & derived()
Definition: EigenBase.h:46
Represents a sequence of transpositions (row/column interchange)
Definition: Transpositions.h:156
static const Eigen::internal::all_t all
Definition: IndexedViewHelper.h:171
ComputationInfo
Definition: Constants.h:440
@ Lower
Definition: Constants.h:209
@ Upper
Definition: Constants.h:211
@ NumericalIssue
Definition: Constants.h:444
@ Success
Definition: Constants.h:442
Namespace containing all symbols from the Eigen library.
Definition: Core:141
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:74
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
Definition: EigenBase.h:30
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:63
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:39
Derived & derived()
Definition: EigenBase.h:46
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:60
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:67