19 template<
typename MatrixType,
int UpLo>
struct LDLT_Traits;
22 enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
50 template<
typename _MatrixType,
int _UpLo>
class LDLT
53 typedef _MatrixType MatrixType;
55 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
58 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
61 typedef typename MatrixType::Scalar Scalar;
64 typedef typename MatrixType::StorageIndex StorageIndex;
70 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
80 m_sign(internal::ZeroSign),
81 m_isInitialized(false)
91 : m_matrix(size, size),
92 m_transpositions(size),
94 m_sign(internal::ZeroSign),
95 m_isInitialized(false)
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)
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)
137 m_isInitialized =
false;
141 inline typename Traits::MatrixU
matrixU()
const
143 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
144 return Traits::getU(m_matrix);
148 inline typename Traits::MatrixL
matrixL()
const
150 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
151 return Traits::getL(m_matrix);
158 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
159 return m_transpositions;
165 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
166 return m_matrix.diagonal();
172 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
173 return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
179 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
180 return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
198 template<
typename Rhs>
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");
208 template<
typename Derived>
211 template<
typename InputType>
219 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
220 return internal::rcond_estimate_helper(m_l1_norm, *
this);
223 template <
typename Derived>
232 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
245 inline Index rows()
const {
return m_matrix.rows(); }
246 inline Index cols()
const {
return m_matrix.cols(); }
255 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
259 #ifndef EIGEN_PARSED_BY_DOXYGEN
260 template<
typename RhsType,
typename DstType>
262 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
267 static void check_template_parameters()
269 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
279 RealScalar m_l1_norm;
280 TranspositionType m_transpositions;
281 TmpMatrixType m_temporary;
282 internal::SignMatrix m_sign;
283 bool m_isInitialized;
289 template<
int UpLo>
struct ldlt_inplace;
291 template<>
struct ldlt_inplace<
Lower>
293 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
294 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix&
sign)
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;
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;
315 for (
Index k = 0; k < size; ++k)
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;
322 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
323 if(k != index_of_biggest_in_corner)
327 Index s = size-index_of_biggest_in_corner-1;
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)
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);
337 if(NumTraits<Scalar>::IsComplex)
338 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
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);
352 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
353 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
355 A21.noalias() -= A20 * temp.head(k);
362 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
363 bool pivot_is_valid = (
abs(realAkk) > RealScalar(0));
365 if(k==0 && !pivot_is_valid)
370 for(
Index j = 0; j<size; ++j)
372 transpositions.coeffRef(j) = IndexType(j);
373 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
378 if((rs>0) && pivot_is_valid)
381 ret = ret && (A21.array()==Scalar(0)).all();
383 if(found_zero_pivot && pivot_is_valid) ret =
false;
384 else if(!pivot_is_valid) found_zero_pivot =
true;
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;
406 template<
typename MatrixType,
typename WDerived>
407 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w,
const typename MatrixType::RealScalar& sigma=1)
409 using numext::isfinite;
410 typedef typename MatrixType::Scalar Scalar;
411 typedef typename MatrixType::RealScalar RealScalar;
413 const Index size = mat.rows();
414 eigen_assert(mat.cols() == size && w.size()==size);
416 RealScalar alpha = 1;
419 for (
Index j = 0; j < size; j++)
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;
431 mat.coeffRef(j,j) += swj2/alpha;
437 w.tail(rs) -= wj * mat.col(j).tail(rs);
439 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
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)
448 tmp = transpositions * w;
450 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
454 template<>
struct ldlt_inplace<
Upper>
456 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
457 static EIGEN_STRONG_INLINE
bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix&
sign)
459 Transpose<MatrixType> matt(mat);
460 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp,
sign);
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)
466 Transpose<MatrixType> matt(mat);
467 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
471 template<
typename MatrixType>
struct LDLT_Traits<MatrixType,
Lower>
473 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
474 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
475 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m); }
476 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m.adjoint()); }
479 template<
typename MatrixType>
struct LDLT_Traits<MatrixType,
Upper>
481 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
482 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
483 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m.adjoint()); }
484 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m); }
491 template<
typename MatrixType,
int _UpLo>
492 template<
typename InputType>
495 check_template_parameters();
503 m_l1_norm = RealScalar(0);
505 for (
Index col = 0; col < size; ++col) {
506 RealScalar abs_col_sum;
508 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
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;
515 m_transpositions.resize(size);
516 m_isInitialized =
false;
517 m_temporary.resize(size);
518 m_sign = internal::ZeroSign;
520 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ?
Success :
NumericalIssue;
522 m_isInitialized =
true;
531 template<
typename MatrixType,
int _UpLo>
532 template<
typename Derived>
535 typedef typename TranspositionType::StorageIndex IndexType;
539 eigen_assert(m_matrix.rows()==size);
543 m_matrix.resize(size,size);
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);
549 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
550 m_isInitialized =
true;
553 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
558 #ifndef EIGEN_PARSED_BY_DOXYGEN
559 template<
typename _MatrixType,
int _UpLo>
560 template<
typename RhsType,
typename DstType>
563 eigen_assert(rhs.rows() == rows());
565 dst = m_transpositions * rhs;
568 matrixL().solveInPlace(dst);
581 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
583 for (
Index i = 0; i < vecD.size(); ++i)
585 if(
abs(vecD(i)) > tolerance)
586 dst.row(i) /= vecD(i);
588 dst.row(i).setZero();
592 matrixU().solveInPlace(dst);
595 dst = m_transpositions.transpose() * dst;
612 template<
typename MatrixType,
int _UpLo>
613 template<
typename Derived>
614 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX)
const
616 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
617 eigen_assert(m_matrix.rows() == bAndX.rows());
619 bAndX = this->solve(bAndX);
627 template<
typename MatrixType,
int _UpLo>
630 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
631 const Index size = m_matrix.rows();
632 MatrixType res(size,size);
636 res = transpositionsP() * res;
638 res = matrixU() * res;
640 res = vectorD().real().asDiagonal() * res;
642 res = matrixL() * res;
644 res = transpositionsP().transpose() * res;
653 template<
typename MatrixType,
unsigned int UpLo>
664 template<
typename Derived>
Derived & derived()
Definition: EigenBase.h:45
Index rows() const
Definition: EigenBase.h:59
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:65
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:51
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
LDLT()
Default Constructor.
Definition: LDLT.h:77
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
bool isPositive() const
Definition: LDLT.h:170
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
void setZero()
Definition: LDLT.h:135
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
bool isNegative(void) const
Definition: LDLT.h:177
const LDLT & adjoint() const
Definition: LDLT.h:243
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
Eigen::Index Index
Definition: LDLT.h:63
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
MatrixType reconstructedMatrix() const
Definition: LDLT.h:628
RealScalar rcond() const
Definition: LDLT.h:217
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:666
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:655
Pseudo expression representing a solving operation.
Definition: Solve.h:63
ComputationInfo
Definition: Constants.h:430
@ Lower
Definition: Constants.h:204
@ Upper
Definition: Constants.h:206
@ NumericalIssue
Definition: Constants.h:434
@ Success
Definition: Constants.h:432
Namespace containing all symbols from the Eigen library.
Definition: Core:309
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
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_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)
Definition: EigenBase.h:30
Index cols() const
Definition: EigenBase.h:62
Derived & derived()
Definition: EigenBase.h:45
Index rows() const
Definition: EigenBase.h:59
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151