Warning, file /include/eigen3/Eigen/src/Cholesky/LDLT.h was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #ifndef EIGEN_LDLT_H
0014 #define EIGEN_LDLT_H
0015
0016 namespace Eigen {
0017
0018 namespace internal {
0019 template<typename _MatrixType, int _UpLo> struct traits<LDLT<_MatrixType, _UpLo> >
0020 : traits<_MatrixType>
0021 {
0022 typedef MatrixXpr XprKind;
0023 typedef SolverStorage StorageKind;
0024 typedef int StorageIndex;
0025 enum { Flags = 0 };
0026 };
0027
0028 template<typename MatrixType, int UpLo> struct LDLT_Traits;
0029
0030
0031 enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
0032 }
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059 template<typename _MatrixType, int _UpLo> class LDLT
0060 : public SolverBase<LDLT<_MatrixType, _UpLo> >
0061 {
0062 public:
0063 typedef _MatrixType MatrixType;
0064 typedef SolverBase<LDLT> Base;
0065 friend class SolverBase<LDLT>;
0066
0067 EIGEN_GENERIC_PUBLIC_INTERFACE(LDLT)
0068 enum {
0069 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0070 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
0071 UpLo = _UpLo
0072 };
0073 typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> TmpMatrixType;
0074
0075 typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
0076 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
0077
0078 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
0079
0080
0081
0082
0083
0084
0085 LDLT()
0086 : m_matrix(),
0087 m_transpositions(),
0088 m_sign(internal::ZeroSign),
0089 m_isInitialized(false)
0090 {}
0091
0092
0093
0094
0095
0096
0097
0098 explicit LDLT(Index size)
0099 : m_matrix(size, size),
0100 m_transpositions(size),
0101 m_temporary(size),
0102 m_sign(internal::ZeroSign),
0103 m_isInitialized(false)
0104 {}
0105
0106
0107
0108
0109
0110
0111
0112 template<typename InputType>
0113 explicit LDLT(const EigenBase<InputType>& matrix)
0114 : m_matrix(matrix.rows(), matrix.cols()),
0115 m_transpositions(matrix.rows()),
0116 m_temporary(matrix.rows()),
0117 m_sign(internal::ZeroSign),
0118 m_isInitialized(false)
0119 {
0120 compute(matrix.derived());
0121 }
0122
0123
0124
0125
0126
0127
0128
0129 template<typename InputType>
0130 explicit LDLT(EigenBase<InputType>& matrix)
0131 : m_matrix(matrix.derived()),
0132 m_transpositions(matrix.rows()),
0133 m_temporary(matrix.rows()),
0134 m_sign(internal::ZeroSign),
0135 m_isInitialized(false)
0136 {
0137 compute(matrix.derived());
0138 }
0139
0140
0141
0142
0143 void setZero()
0144 {
0145 m_isInitialized = false;
0146 }
0147
0148
0149 inline typename Traits::MatrixU matrixU() const
0150 {
0151 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0152 return Traits::getU(m_matrix);
0153 }
0154
0155
0156 inline typename Traits::MatrixL matrixL() const
0157 {
0158 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0159 return Traits::getL(m_matrix);
0160 }
0161
0162
0163
0164 inline const TranspositionType& transpositionsP() const
0165 {
0166 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0167 return m_transpositions;
0168 }
0169
0170
0171 inline Diagonal<const MatrixType> vectorD() const
0172 {
0173 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0174 return m_matrix.diagonal();
0175 }
0176
0177
0178 inline bool isPositive() const
0179 {
0180 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0181 return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
0182 }
0183
0184
0185 inline bool isNegative(void) const
0186 {
0187 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0188 return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
0189 }
0190
0191 #ifdef EIGEN_PARSED_BY_DOXYGEN
0192
0193
0194
0195
0196
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206
0207 template<typename Rhs>
0208 inline const Solve<LDLT, Rhs>
0209 solve(const MatrixBase<Rhs>& b) const;
0210 #endif
0211
0212 template<typename Derived>
0213 bool solveInPlace(MatrixBase<Derived> &bAndX) const;
0214
0215 template<typename InputType>
0216 LDLT& compute(const EigenBase<InputType>& matrix);
0217
0218
0219
0220
0221 RealScalar rcond() const
0222 {
0223 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0224 return internal::rcond_estimate_helper(m_l1_norm, *this);
0225 }
0226
0227 template <typename Derived>
0228 LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
0229
0230
0231
0232
0233
0234 inline const MatrixType& matrixLDLT() const
0235 {
0236 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0237 return m_matrix;
0238 }
0239
0240 MatrixType reconstructedMatrix() const;
0241
0242
0243
0244
0245
0246
0247 const LDLT& adjoint() const { return *this; };
0248
0249 EIGEN_DEVICE_FUNC inline EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); }
0250 EIGEN_DEVICE_FUNC inline EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
0251
0252
0253
0254
0255
0256
0257 ComputationInfo info() const
0258 {
0259 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0260 return m_info;
0261 }
0262
0263 #ifndef EIGEN_PARSED_BY_DOXYGEN
0264 template<typename RhsType, typename DstType>
0265 void _solve_impl(const RhsType &rhs, DstType &dst) const;
0266
0267 template<bool Conjugate, typename RhsType, typename DstType>
0268 void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
0269 #endif
0270
0271 protected:
0272
0273 static void check_template_parameters()
0274 {
0275 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
0276 }
0277
0278
0279
0280
0281
0282
0283
0284 MatrixType m_matrix;
0285 RealScalar m_l1_norm;
0286 TranspositionType m_transpositions;
0287 TmpMatrixType m_temporary;
0288 internal::SignMatrix m_sign;
0289 bool m_isInitialized;
0290 ComputationInfo m_info;
0291 };
0292
0293 namespace internal {
0294
0295 template<int UpLo> struct ldlt_inplace;
0296
0297 template<> struct ldlt_inplace<Lower>
0298 {
0299 template<typename MatrixType, typename TranspositionType, typename Workspace>
0300 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
0301 {
0302 using std::abs;
0303 typedef typename MatrixType::Scalar Scalar;
0304 typedef typename MatrixType::RealScalar RealScalar;
0305 typedef typename TranspositionType::StorageIndex IndexType;
0306 eigen_assert(mat.rows()==mat.cols());
0307 const Index size = mat.rows();
0308 bool found_zero_pivot = false;
0309 bool ret = true;
0310
0311 if (size <= 1)
0312 {
0313 transpositions.setIdentity();
0314 if(size==0) sign = ZeroSign;
0315 else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
0316 else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
0317 else sign = ZeroSign;
0318 return true;
0319 }
0320
0321 for (Index k = 0; k < size; ++k)
0322 {
0323
0324 Index index_of_biggest_in_corner;
0325 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
0326 index_of_biggest_in_corner += k;
0327
0328 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
0329 if(k != index_of_biggest_in_corner)
0330 {
0331
0332
0333 Index s = size-index_of_biggest_in_corner-1;
0334 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
0335 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
0336 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
0337 for(Index i=k+1;i<index_of_biggest_in_corner;++i)
0338 {
0339 Scalar tmp = mat.coeffRef(i,k);
0340 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
0341 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
0342 }
0343 if(NumTraits<Scalar>::IsComplex)
0344 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
0345 }
0346
0347
0348
0349
0350
0351 Index rs = size - k - 1;
0352 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
0353 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
0354 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
0355
0356 if(k>0)
0357 {
0358 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
0359 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
0360 if(rs>0)
0361 A21.noalias() -= A20 * temp.head(k);
0362 }
0363
0364
0365
0366
0367
0368 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
0369 bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
0370
0371 if(k==0 && !pivot_is_valid)
0372 {
0373
0374
0375 sign = ZeroSign;
0376 for(Index j = 0; j<size; ++j)
0377 {
0378 transpositions.coeffRef(j) = IndexType(j);
0379 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
0380 }
0381 return ret;
0382 }
0383
0384 if((rs>0) && pivot_is_valid)
0385 A21 /= realAkk;
0386 else if(rs>0)
0387 ret = ret && (A21.array()==Scalar(0)).all();
0388
0389 if(found_zero_pivot && pivot_is_valid) ret = false;
0390 else if(!pivot_is_valid) found_zero_pivot = true;
0391
0392 if (sign == PositiveSemiDef) {
0393 if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
0394 } else if (sign == NegativeSemiDef) {
0395 if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
0396 } else if (sign == ZeroSign) {
0397 if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
0398 else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
0399 }
0400 }
0401
0402 return ret;
0403 }
0404
0405
0406
0407
0408
0409
0410
0411
0412 template<typename MatrixType, typename WDerived>
0413 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
0414 {
0415 using numext::isfinite;
0416 typedef typename MatrixType::Scalar Scalar;
0417 typedef typename MatrixType::RealScalar RealScalar;
0418
0419 const Index size = mat.rows();
0420 eigen_assert(mat.cols() == size && w.size()==size);
0421
0422 RealScalar alpha = 1;
0423
0424
0425 for (Index j = 0; j < size; j++)
0426 {
0427
0428 if (!(isfinite)(alpha))
0429 break;
0430
0431
0432 RealScalar dj = numext::real(mat.coeff(j,j));
0433 Scalar wj = w.coeff(j);
0434 RealScalar swj2 = sigma*numext::abs2(wj);
0435 RealScalar gamma = dj*alpha + swj2;
0436
0437 mat.coeffRef(j,j) += swj2/alpha;
0438 alpha += swj2/dj;
0439
0440
0441
0442 Index rs = size-j-1;
0443 w.tail(rs) -= wj * mat.col(j).tail(rs);
0444 if(gamma != 0)
0445 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
0446 }
0447 return true;
0448 }
0449
0450 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
0451 static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
0452 {
0453
0454 tmp = transpositions * w;
0455
0456 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
0457 }
0458 };
0459
0460 template<> struct ldlt_inplace<Upper>
0461 {
0462 template<typename MatrixType, typename TranspositionType, typename Workspace>
0463 static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
0464 {
0465 Transpose<MatrixType> matt(mat);
0466 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
0467 }
0468
0469 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
0470 static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
0471 {
0472 Transpose<MatrixType> matt(mat);
0473 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
0474 }
0475 };
0476
0477 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
0478 {
0479 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
0480 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
0481 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
0482 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
0483 };
0484
0485 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
0486 {
0487 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
0488 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
0489 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
0490 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
0491 };
0492
0493 }
0494
0495
0496
0497 template<typename MatrixType, int _UpLo>
0498 template<typename InputType>
0499 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
0500 {
0501 check_template_parameters();
0502
0503 eigen_assert(a.rows()==a.cols());
0504 const Index size = a.rows();
0505
0506 m_matrix = a.derived();
0507
0508
0509 m_l1_norm = RealScalar(0);
0510
0511 for (Index col = 0; col < size; ++col) {
0512 RealScalar abs_col_sum;
0513 if (_UpLo == Lower)
0514 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
0515 else
0516 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
0517 if (abs_col_sum > m_l1_norm)
0518 m_l1_norm = abs_col_sum;
0519 }
0520
0521 m_transpositions.resize(size);
0522 m_isInitialized = false;
0523 m_temporary.resize(size);
0524 m_sign = internal::ZeroSign;
0525
0526 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
0527
0528 m_isInitialized = true;
0529 return *this;
0530 }
0531
0532
0533
0534
0535
0536
0537 template<typename MatrixType, int _UpLo>
0538 template<typename Derived>
0539 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
0540 {
0541 typedef typename TranspositionType::StorageIndex IndexType;
0542 const Index size = w.rows();
0543 if (m_isInitialized)
0544 {
0545 eigen_assert(m_matrix.rows()==size);
0546 }
0547 else
0548 {
0549 m_matrix.resize(size,size);
0550 m_matrix.setZero();
0551 m_transpositions.resize(size);
0552 for (Index i = 0; i < size; i++)
0553 m_transpositions.coeffRef(i) = IndexType(i);
0554 m_temporary.resize(size);
0555 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
0556 m_isInitialized = true;
0557 }
0558
0559 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
0560
0561 return *this;
0562 }
0563
0564 #ifndef EIGEN_PARSED_BY_DOXYGEN
0565 template<typename _MatrixType, int _UpLo>
0566 template<typename RhsType, typename DstType>
0567 void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
0568 {
0569 _solve_impl_transposed<true>(rhs, dst);
0570 }
0571
0572 template<typename _MatrixType,int _UpLo>
0573 template<bool Conjugate, typename RhsType, typename DstType>
0574 void LDLT<_MatrixType,_UpLo>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
0575 {
0576
0577 dst = m_transpositions * rhs;
0578
0579
0580
0581 matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
0582
0583
0584
0585
0586 using std::abs;
0587 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
0588
0589
0590
0591
0592
0593
0594
0595 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
0596 for (Index i = 0; i < vecD.size(); ++i)
0597 {
0598 if(abs(vecD(i)) > tolerance)
0599 dst.row(i) /= vecD(i);
0600 else
0601 dst.row(i).setZero();
0602 }
0603
0604
0605
0606 matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
0607
0608
0609
0610 dst = m_transpositions.transpose() * dst;
0611 }
0612 #endif
0613
0614
0615
0616
0617
0618
0619
0620
0621
0622
0623
0624
0625
0626
0627 template<typename MatrixType,int _UpLo>
0628 template<typename Derived>
0629 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
0630 {
0631 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0632 eigen_assert(m_matrix.rows() == bAndX.rows());
0633
0634 bAndX = this->solve(bAndX);
0635
0636 return true;
0637 }
0638
0639
0640
0641
0642 template<typename MatrixType, int _UpLo>
0643 MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
0644 {
0645 eigen_assert(m_isInitialized && "LDLT is not initialized.");
0646 const Index size = m_matrix.rows();
0647 MatrixType res(size,size);
0648
0649
0650 res.setIdentity();
0651 res = transpositionsP() * res;
0652
0653 res = matrixU() * res;
0654
0655 res = vectorD().real().asDiagonal() * res;
0656
0657 res = matrixL() * res;
0658
0659 res = transpositionsP().transpose() * res;
0660
0661 return res;
0662 }
0663
0664
0665
0666
0667
0668 template<typename MatrixType, unsigned int UpLo>
0669 inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
0670 SelfAdjointView<MatrixType, UpLo>::ldlt() const
0671 {
0672 return LDLT<PlainObject,UpLo>(m_matrix);
0673 }
0674
0675
0676
0677
0678
0679 template<typename Derived>
0680 inline const LDLT<typename MatrixBase<Derived>::PlainObject>
0681 MatrixBase<Derived>::ldlt() const
0682 {
0683 return LDLT<PlainObject>(derived());
0684 }
0685
0686 }
0687
0688 #endif