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;
81 m_isInitialized(false)
91 : m_matrix(size, size),
92 m_transpositions(size),
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()),
110 m_isInitialized(false)
121 template<
typename InputType>
123 : m_matrix(matrix.derived()),
124 m_transpositions(matrix.rows()),
125 m_temporary(matrix.rows()),
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.");
223 template <
typename Derived>
232 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
236 MatrixType reconstructedMatrix()
const;
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;
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 (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
309 else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
310 else sign = ZeroSign;
314 for (
Index k = 0; k < size; ++k)
317 Index index_of_biggest_in_corner;
318 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
319 index_of_biggest_in_corner += k;
321 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
322 if(k != index_of_biggest_in_corner)
326 Index s = size-index_of_biggest_in_corner-1;
327 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
328 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
329 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
330 for(
Index i=k+1;i<index_of_biggest_in_corner;++i)
332 Scalar tmp = mat.coeffRef(i,k);
333 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
334 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
337 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
344 Index rs = size - k - 1;
351 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
352 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
354 A21.noalias() -= A20 * temp.head(k);
361 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
362 bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
364 if(k==0 && !pivot_is_valid)
369 for(
Index j = 0; j<size; ++j)
371 transpositions.coeffRef(j) = IndexType(j);
372 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
377 if((rs>0) && pivot_is_valid)
380 if(found_zero_pivot && pivot_is_valid) ret =
false;
381 else if(!pivot_is_valid) found_zero_pivot =
true;
383 if (sign == PositiveSemiDef) {
384 if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
385 }
else if (sign == NegativeSemiDef) {
386 if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
387 }
else if (sign == ZeroSign) {
388 if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
389 else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
403 template<
typename MatrixType,
typename WDerived>
404 static bool updateInPlace(MatrixType& mat,
MatrixBase<WDerived>& w,
const typename MatrixType::RealScalar& sigma=1)
406 using numext::isfinite;
407 typedef typename MatrixType::Scalar Scalar;
408 typedef typename MatrixType::RealScalar RealScalar;
410 const Index size = mat.rows();
411 eigen_assert(mat.cols() == size && w.size()==size);
413 RealScalar alpha = 1;
416 for (
Index j = 0; j < size; j++)
419 if (!(isfinite)(alpha))
423 RealScalar dj = numext::real(mat.coeff(j,j));
424 Scalar wj = w.coeff(j);
425 RealScalar swj2 = sigma*numext::abs2(wj);
426 RealScalar gamma = dj*alpha + swj2;
428 mat.coeffRef(j,j) += swj2/alpha;
434 w.
tail(rs) -= wj * mat.col(j).tail(rs);
436 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.
tail(rs);
441 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
442 static bool update(MatrixType& mat,
const TranspositionType& transpositions, Workspace& tmp,
const WType& w,
const typename MatrixType::RealScalar& sigma=1)
445 tmp = transpositions * w;
453 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
454 static EIGEN_STRONG_INLINE
bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
460 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
461 static EIGEN_STRONG_INLINE
bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w,
const typename MatrixType::RealScalar& sigma=1)
472 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m); }
473 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m.adjoint()); }
480 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m.adjoint()); }
481 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m); }
488 template<
typename MatrixType,
int _UpLo>
489 template<
typename InputType>
492 check_template_parameters();
500 m_l1_norm = RealScalar(0);
502 for (
Index col = 0; col < size; ++col) {
503 RealScalar abs_col_sum;
505 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
507 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
508 if (abs_col_sum > m_l1_norm)
509 m_l1_norm = abs_col_sum;
512 m_transpositions.resize(size);
513 m_isInitialized =
false;
514 m_temporary.resize(size);
515 m_sign = internal::ZeroSign;
519 m_isInitialized =
true;
528 template<
typename MatrixType,
int _UpLo>
529 template<
typename Derived>
532 typedef typename TranspositionType::StorageIndex IndexType;
533 const Index size = w.rows();
536 eigen_assert(m_matrix.rows()==size);
540 m_matrix.resize(size,size);
542 m_transpositions.resize(size);
543 for (
Index i = 0; i < size; i++)
544 m_transpositions.coeffRef(i) = IndexType(i);
545 m_temporary.resize(size);
546 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
547 m_isInitialized =
true;
555 #ifndef EIGEN_PARSED_BY_DOXYGEN 556 template<
typename _MatrixType,
int _UpLo>
557 template<
typename RhsType,
typename DstType>
560 eigen_assert(rhs.rows() == rows());
562 dst = m_transpositions * rhs;
565 matrixL().solveInPlace(dst);
579 for (
Index i = 0; i < vecD.size(); ++i)
581 if(abs(vecD(i)) > tolerance)
582 dst.row(i) /= vecD(i);
584 dst.row(i).setZero();
588 matrixU().solveInPlace(dst);
591 dst = m_transpositions.transpose() * dst;
608 template<
typename MatrixType,
int _UpLo>
609 template<
typename Derived>
612 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
613 eigen_assert(m_matrix.rows() == bAndX.rows());
615 bAndX = this->solve(bAndX);
623 template<
typename MatrixType,
int _UpLo>
626 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
627 const Index size = m_matrix.rows();
628 MatrixType res(size,size);
632 res = transpositionsP() * res;
634 res = matrixU() * res;
636 res = vectorD().real().asDiagonal() * res;
638 res = matrixL() * res;
640 res = transpositionsP().transpose() * res;
649 template<
typename MatrixType,
unsigned int UpLo>
660 template<
typename Derived>
669 #endif // EIGEN_LDLT_H ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:50
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:651
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
EIGEN_DEVICE_FUNC Index rows() const
Definition: EigenBase.h:58
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
RealScalar rcond() const
Definition: LDLT.h:217
Expression of the transpose of a matrix.
Definition: Transpose.h:52
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
Definition: ConditionEstimator.h:159
void setZero()
Clear any existing decomposition.
Definition: LDLT.h:135
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
bool isNegative(void) const
Definition: LDLT.h:177
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
LDLT()
Default Constructor.
Definition: LDLT.h:77
View matrix as a lower triangular matrix.
Definition: Constants.h:204
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
The provided data did not satisfy the prerequisites.
Definition: Constants.h:434
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
const LDLT & adjoint() const
Definition: LDLT.h:243
View matrix as an upper triangular matrix.
Definition: Constants.h:206
Computation was successful.
Definition: Constants.h:432
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:662
bool isPositive() const
Definition: LDLT.h:170
Eigen::Index Index
Definition: LDLT.h:63
Definition: BandTriangularSolver.h:13
EIGEN_DEVICE_FUNC SegmentReturnType tail(Index n)
Definition: DenseBase.h:950
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:186
EIGEN_DEVICE_FUNC Index cols() const
Definition: EigenBase.h:61
MatrixType reconstructedMatrix() const
Definition: LDLT.h:624
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
Pseudo expression representing a solving operation.
Definition: Solve.h:62
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:44
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48