13#include "./InternalHeaderCheck.h"
19template<
typename MatrixType_,
int UpLo_>
struct traits<LLT<MatrixType_, UpLo_> >
22 typedef MatrixXpr XprKind;
23 typedef SolverStorage StorageKind;
24 typedef int StorageIndex;
28template<
typename MatrixType,
int UpLo>
struct LLT_Traits;
68template<
typename MatrixType_,
int UpLo_>
class LLT
72 typedef MatrixType_ MatrixType;
76 EIGEN_GENERIC_PUBLIC_INTERFACE(
LLT)
78 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
82 PacketSize = internal::packet_traits<Scalar>::size,
83 AlignmentMask = int(PacketSize)-1,
87 typedef internal::LLT_Traits<MatrixType,UpLo> Traits;
95 LLT() : m_matrix(), m_isInitialized(false) {}
104 m_isInitialized(false) {}
106 template<
typename InputType>
108 : m_matrix(matrix.rows(), matrix.cols()),
109 m_isInitialized(false)
121 template<
typename InputType>
124 m_isInitialized(false)
130 inline typename Traits::MatrixU
matrixU()
const
132 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
133 return Traits::getU(m_matrix);
137 inline typename Traits::MatrixL
matrixL()
const
139 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
140 return Traits::getL(m_matrix);
143 #ifdef EIGEN_PARSED_BY_DOXYGEN
154 template<
typename Rhs>
159 template<
typename Derived>
162 template<
typename InputType>
170 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
171 eigen_assert(m_info ==
Success &&
"LLT failed because matrix appears to be negative");
172 return internal::rcond_estimate_helper(m_l1_norm, *
this);
181 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
195 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
206 inline EIGEN_CONSTEXPR
Index rows() const EIGEN_NOEXCEPT {
return m_matrix.rows(); }
207 inline EIGEN_CONSTEXPR
Index cols() const EIGEN_NOEXCEPT {
return m_matrix.cols(); }
209 template<
typename VectorType>
210 LLT & rankUpdate(
const VectorType& vec,
const RealScalar& sigma = 1);
212 #ifndef EIGEN_PARSED_BY_DOXYGEN
213 template<
typename RhsType,
typename DstType>
214 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
216 template<
bool Conjugate,
typename RhsType,
typename DstType>
217 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
222 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
229 RealScalar m_l1_norm;
230 bool m_isInitialized;
236template<
typename Scalar,
int UpLo>
struct llt_inplace;
238template<
typename MatrixType,
typename VectorType>
239static Index llt_rank_update_lower(MatrixType& mat,
const VectorType& vec,
const typename MatrixType::RealScalar& sigma)
242 typedef typename MatrixType::Scalar Scalar;
243 typedef typename MatrixType::RealScalar RealScalar;
244 typedef typename MatrixType::ColXpr ColXpr;
245 typedef typename internal::remove_all<ColXpr>::type ColXprCleaned;
246 typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
247 typedef Matrix<Scalar,Dynamic,1> TempVectorType;
248 typedef typename TempVectorType::SegmentReturnType TempVecSegment;
250 Index n = mat.cols();
251 eigen_assert(mat.rows()==n && vec.size()==n);
260 temp =
sqrt(sigma) * vec;
262 for(
Index i=0; i<n; ++i)
264 JacobiRotation<Scalar> g;
265 g.makeGivens(mat(i,i), -temp(i), &mat(i,i));
270 ColXprSegment x(mat.col(i).tail(rs));
271 TempVecSegment y(temp.tail(rs));
272 apply_rotation_in_the_plane(x, y, g);
280 for(
Index j=0; j<n; ++j)
282 RealScalar Ljj = numext::real(mat.coeff(j,j));
283 RealScalar dj = numext::abs2(Ljj);
284 Scalar wj = temp.coeff(j);
285 RealScalar swj2 = sigma*numext::abs2(wj);
286 RealScalar gamma = dj*beta + swj2;
288 RealScalar x = dj + swj2/beta;
289 if (x<=RealScalar(0))
291 RealScalar nLjj =
sqrt(x);
292 mat.coeffRef(j,j) = nLjj;
299 temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
301 mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
308template<
typename Scalar>
struct llt_inplace<Scalar,
Lower>
310 typedef typename NumTraits<Scalar>::Real RealScalar;
311 template<
typename MatrixType>
312 static Index unblocked(MatrixType& mat)
316 eigen_assert(mat.rows()==mat.cols());
317 const Index size = mat.rows();
318 for(
Index k = 0; k < size; ++k)
322 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
323 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
324 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
326 RealScalar x = numext::real(mat.coeff(k,k));
327 if (k>0) x -= A10.squaredNorm();
328 if (x<=RealScalar(0))
330 mat.coeffRef(k,k) = x =
sqrt(x);
331 if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
337 template<
typename MatrixType>
338 static Index blocked(MatrixType& m)
340 eigen_assert(m.rows()==m.cols());
341 Index size = m.rows();
345 Index blockSize = size/8;
346 blockSize = (blockSize/16)*16;
347 blockSize = (std::min)((std::max)(blockSize,
Index(8)),
Index(128));
349 for (
Index k=0; k<size; k+=blockSize)
355 Index bs = (std::min)(blockSize, size-k);
362 if((ret=unblocked(A11))>=0)
return k+ret;
363 if(rs>0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
369 template<
typename MatrixType,
typename VectorType>
370 static Index rankUpdate(MatrixType& mat,
const VectorType& vec,
const RealScalar& sigma)
372 return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
376template<
typename Scalar>
struct llt_inplace<Scalar,
Upper>
378 typedef typename NumTraits<Scalar>::Real RealScalar;
380 template<
typename MatrixType>
381 static EIGEN_STRONG_INLINE
Index unblocked(MatrixType& mat)
383 Transpose<MatrixType> matt(mat);
384 return llt_inplace<Scalar, Lower>::unblocked(matt);
386 template<
typename MatrixType>
387 static EIGEN_STRONG_INLINE
Index blocked(MatrixType& mat)
389 Transpose<MatrixType> matt(mat);
390 return llt_inplace<Scalar, Lower>::blocked(matt);
392 template<
typename MatrixType,
typename VectorType>
393 static Index rankUpdate(MatrixType& mat,
const VectorType& vec,
const RealScalar& sigma)
395 Transpose<MatrixType> matt(mat);
396 return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
400template<
typename MatrixType>
struct LLT_Traits<MatrixType,
Lower>
402 typedef const TriangularView<const MatrixType, Lower> MatrixL;
403 typedef const TriangularView<const typename MatrixType::AdjointReturnType, Upper> MatrixU;
404 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m); }
405 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m.adjoint()); }
406 static bool inplace_decomposition(MatrixType& m)
407 {
return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m)==-1; }
410template<
typename MatrixType>
struct LLT_Traits<MatrixType,
Upper>
412 typedef const TriangularView<const typename MatrixType::AdjointReturnType, Lower> MatrixL;
413 typedef const TriangularView<const MatrixType, Upper> MatrixU;
414 static inline MatrixL getL(
const MatrixType& m) {
return MatrixL(m.adjoint()); }
415 static inline MatrixU getU(
const MatrixType& m) {
return MatrixU(m); }
416 static bool inplace_decomposition(MatrixType& m)
417 {
return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m)==-1; }
429template<
typename MatrixType,
int UpLo_>
430template<
typename InputType>
435 m_matrix.resize(size, size);
436 if (!internal::is_same_dense(m_matrix, a.
derived()))
440 m_l1_norm = RealScalar(0);
442 for (
Index col = 0; col < size; ++col) {
443 RealScalar abs_col_sum;
445 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
447 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
448 if (abs_col_sum > m_l1_norm)
449 m_l1_norm = abs_col_sum;
452 m_isInitialized =
true;
453 bool ok = Traits::inplace_decomposition(m_matrix);
464template<
typename MatrixType_,
int UpLo_>
465template<
typename VectorType>
468 EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
469 eigen_assert(v.size()==m_matrix.cols());
470 eigen_assert(m_isInitialized);
471 if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
479#ifndef EIGEN_PARSED_BY_DOXYGEN
480template<
typename MatrixType_,
int UpLo_>
481template<
typename RhsType,
typename DstType>
484 _solve_impl_transposed<true>(rhs, dst);
487template<
typename MatrixType_,
int UpLo_>
488template<
bool Conjugate,
typename RhsType,
typename DstType>
489void LLT<MatrixType_,UpLo_>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const
493 matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
494 matrixU().template conjugateIf<!Conjugate>().solveInPlace(dst);
511template<
typename MatrixType,
int UpLo_>
512template<
typename Derived>
513void LLT<MatrixType,UpLo_>::solveInPlace(
const MatrixBase<Derived> &bAndX)
const
515 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
516 eigen_assert(m_matrix.rows()==bAndX.rows());
517 matrixL().solveInPlace(bAndX);
518 matrixU().solveInPlace(bAndX);
524template<
typename MatrixType,
int UpLo_>
527 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
528 return matrixL() * matrixL().adjoint().toDenseMatrix();
535template<
typename Derived>
546template<
typename MatrixType,
unsigned int UpLo>
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:107
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LLT.h:193
LLT(Index size)
Default Constructor with memory preallocation.
Definition: LLT.h:103
RealScalar rcond() const
Definition: LLT.h:168
const MatrixType & matrixLLT() const
Definition: LLT.h:179
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
LLT(EigenBase< InputType > &matrix)
Constructs a LLT factorization from a given matrix.
Definition: LLT.h:122
Traits::MatrixU matrixU() const
Definition: LLT.h:130
const LLT & adjoint() const EIGEN_NOEXCEPT
Definition: LLT.h:204
LLT()
Default Constructor.
Definition: LLT.h:95
MatrixType reconstructedMatrix() const
Definition: LLT.h:525
Traits::MatrixL matrixL() const
Definition: LLT.h:137
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const LLT< PlainObject > llt() const
Definition: LLT.h:537
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:548
Pseudo expression representing a solving operation.
Definition: Solve.h:65
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:71
LLT< MatrixType_, UpLo_ > & derived()
Definition: EigenBase.h:48
ComputationInfo
Definition: Constants.h:442
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:446
@ Success
Definition: Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
Definition: EigenBase.h:32
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Derived & derived()
Definition: EigenBase.h:48
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:235