LLT.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_LLT_H
11 #define EIGEN_LLT_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal{
18 
19 template<typename MatrixType_, int UpLo_> struct traits<LLT<MatrixType_, UpLo_> >
20  : traits<MatrixType_>
21 {
22  typedef MatrixXpr XprKind;
23  typedef SolverStorage StorageKind;
24  typedef int StorageIndex;
25  enum { Flags = 0 };
26 };
27 
28 template<typename MatrixType, int UpLo> struct LLT_Traits;
29 }
30 
68 template<typename MatrixType_, int UpLo_> class LLT
69  : public SolverBase<LLT<MatrixType_, UpLo_> >
70 {
71  public:
72  typedef MatrixType_ MatrixType;
74  friend class SolverBase<LLT>;
75 
77  enum {
78  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
79  };
80 
81  enum {
84  UpLo = UpLo_
85  };
86 
87  typedef internal::LLT_Traits<MatrixType,UpLo> Traits;
88 
95  LLT() : m_matrix(), m_isInitialized(false) {}
96 
103  explicit LLT(Index size) : m_matrix(size, size),
104  m_isInitialized(false) {}
105 
106  template<typename InputType>
107  explicit LLT(const EigenBase<InputType>& matrix)
108  : m_matrix(matrix.rows(), matrix.cols()),
109  m_isInitialized(false)
110  {
111  compute(matrix.derived());
112  }
113 
121  template<typename InputType>
122  explicit LLT(EigenBase<InputType>& matrix)
123  : m_matrix(matrix.derived()),
124  m_isInitialized(false)
125  {
126  compute(matrix.derived());
127  }
128 
130  inline typename Traits::MatrixU matrixU() const
131  {
132  eigen_assert(m_isInitialized && "LLT is not initialized.");
133  return Traits::getU(m_matrix);
134  }
135 
137  inline typename Traits::MatrixL matrixL() const
138  {
139  eigen_assert(m_isInitialized && "LLT is not initialized.");
140  return Traits::getL(m_matrix);
141  }
142 
143  #ifdef EIGEN_PARSED_BY_DOXYGEN
154  template<typename Rhs>
155  inline const Solve<LLT, Rhs>
156  solve(const MatrixBase<Rhs>& b) const;
157  #endif
158 
159  template<typename Derived>
160  void solveInPlace(const MatrixBase<Derived> &bAndX) const;
161 
162  template<typename InputType>
164 
168  RealScalar rcond() const
169  {
170  eigen_assert(m_isInitialized && "LLT is not initialized.");
171  eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
173  }
174 
179  inline const MatrixType& matrixLLT() const
180  {
181  eigen_assert(m_isInitialized && "LLT is not initialized.");
182  return m_matrix;
183  }
184 
186 
187 
194  {
195  eigen_assert(m_isInitialized && "LLT is not initialized.");
196  return m_info;
197  }
198 
204  const LLT& adjoint() const EIGEN_NOEXCEPT { return *this; }
205 
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(); }
208 
209  template<typename VectorType>
210  LLT & rankUpdate(const VectorType& vec, const RealScalar& sigma = 1);
211 
212  #ifndef EIGEN_PARSED_BY_DOXYGEN
213  template<typename RhsType, typename DstType>
214  void _solve_impl(const RhsType &rhs, DstType &dst) const;
215 
216  template<bool Conjugate, typename RhsType, typename DstType>
217  void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
218  #endif
219 
220  protected:
221 
223 
224 
229  RealScalar m_l1_norm;
232 };
233 
234 namespace internal {
235 
236 template<typename Scalar, int UpLo> struct llt_inplace;
237 
238 template<typename MatrixType, typename VectorType>
239 static Index llt_rank_update_lower(MatrixType& mat, const VectorType& vec, const typename MatrixType::RealScalar& sigma)
240 {
241  using std::sqrt;
242  typedef typename MatrixType::Scalar Scalar;
243  typedef typename MatrixType::RealScalar RealScalar;
244  typedef typename MatrixType::ColXpr ColXpr;
245  typedef internal::remove_all_t<ColXpr> ColXprCleaned;
246  typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
247  typedef Matrix<Scalar,Dynamic,1> TempVectorType;
248  typedef typename TempVectorType::SegmentReturnType TempVecSegment;
249 
250  Index n = mat.cols();
251  eigen_assert(mat.rows()==n && vec.size()==n);
252 
253  TempVectorType temp;
254 
255  if(sigma>0)
256  {
257  // This version is based on Givens rotations.
258  // It is faster than the other one below, but only works for updates,
259  // i.e., for sigma > 0
260  temp = sqrt(sigma) * vec;
261 
262  for(Index i=0; i<n; ++i)
263  {
265  g.makeGivens(mat(i,i), -temp(i), &mat(i,i));
266 
267  Index rs = n-i-1;
268  if(rs>0)
269  {
270  ColXprSegment x(mat.col(i).tail(rs));
271  TempVecSegment y(temp.tail(rs));
273  }
274  }
275  }
276  else
277  {
278  temp = vec;
279  RealScalar beta = 1;
280  for(Index j=0; j<n; ++j)
281  {
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;
287 
288  RealScalar x = dj + swj2/beta;
289  if (x<=RealScalar(0))
290  return j;
291  RealScalar nLjj = sqrt(x);
292  mat.coeffRef(j,j) = nLjj;
293  beta += swj2/dj;
294 
295  // Update the terms of L
296  Index rs = n-j-1;
297  if(rs)
298  {
299  temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
300  if(!numext::is_exactly_zero(gamma))
301  mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
302  }
303  }
304  }
305  return -1;
306 }
307 
308 template<typename Scalar> struct llt_inplace<Scalar, Lower>
309 {
310  typedef typename NumTraits<Scalar>::Real RealScalar;
311  template<typename MatrixType>
312  static Index unblocked(MatrixType& mat)
313  {
314  using std::sqrt;
315 
316  eigen_assert(mat.rows()==mat.cols());
317  const Index size = mat.rows();
318  for(Index k = 0; k < size; ++k)
319  {
320  Index rs = size-k-1; // remaining size
321 
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);
325 
326  RealScalar x = numext::real(mat.coeff(k,k));
327  if (k>0) x -= A10.squaredNorm();
328  if (x<=RealScalar(0))
329  return k;
330  mat.coeffRef(k,k) = x = sqrt(x);
331  if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
332  if (rs>0) A21 /= x;
333  }
334  return -1;
335  }
336 
337  template<typename MatrixType>
338  static Index blocked(MatrixType& m)
339  {
340  eigen_assert(m.rows()==m.cols());
341  Index size = m.rows();
342  if(size<32)
343  return unblocked(m);
344 
345  Index blockSize = size/8;
346  blockSize = (blockSize/16)*16;
347  blockSize = (std::min)((std::max)(blockSize,Index(8)), Index(128));
348 
349  for (Index k=0; k<size; k+=blockSize)
350  {
351  // partition the matrix:
352  // A00 | - | -
353  // lu = A10 | A11 | -
354  // A20 | A21 | A22
355  Index bs = (std::min)(blockSize, size-k);
356  Index rs = size - k - bs;
357  Block<MatrixType,Dynamic,Dynamic> A11(m,k, k, bs,bs);
358  Block<MatrixType,Dynamic,Dynamic> A21(m,k+bs,k, rs,bs);
359  Block<MatrixType,Dynamic,Dynamic> A22(m,k+bs,k+bs,rs,rs);
360 
361  Index ret;
362  if((ret=unblocked(A11))>=0) return k+ret;
363  if(rs>0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
364  if(rs>0) A22.template selfadjointView<Lower>().rankUpdate(A21,typename NumTraits<RealScalar>::Literal(-1)); // bottleneck
365  }
366  return -1;
367  }
368 
369  template<typename MatrixType, typename VectorType>
370  static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
371  {
372  return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
373  }
374 };
375 
376 template<typename Scalar> struct llt_inplace<Scalar, Upper>
377 {
378  typedef typename NumTraits<Scalar>::Real RealScalar;
379 
380  template<typename MatrixType>
381  static EIGEN_STRONG_INLINE Index unblocked(MatrixType& mat)
382  {
383  Transpose<MatrixType> matt(mat);
384  return llt_inplace<Scalar, Lower>::unblocked(matt);
385  }
386  template<typename MatrixType>
387  static EIGEN_STRONG_INLINE Index blocked(MatrixType& mat)
388  {
389  Transpose<MatrixType> matt(mat);
390  return llt_inplace<Scalar, Lower>::blocked(matt);
391  }
392  template<typename MatrixType, typename VectorType>
393  static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
394  {
395  Transpose<MatrixType> matt(mat);
396  return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
397  }
398 };
399 
400 template<typename MatrixType> struct LLT_Traits<MatrixType,Lower>
401 {
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; }
408 };
409 
410 template<typename MatrixType> struct LLT_Traits<MatrixType,Upper>
411 {
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; }
418 };
419 
420 } // end namespace internal
421 
429 template<typename MatrixType, int UpLo_>
430 template<typename InputType>
432 {
433  eigen_assert(a.rows()==a.cols());
434  const Index size = a.rows();
435  m_matrix.resize(size, size);
436  if (!internal::is_same_dense(m_matrix, a.derived()))
437  m_matrix = a.derived();
438 
439  // Compute matrix L1 norm = max abs column sum.
440  m_l1_norm = RealScalar(0);
441  // TODO move this code to SelfAdjointView
442  for (Index col = 0; col < size; ++col) {
443  RealScalar abs_col_sum;
444  if (UpLo_ == Lower)
445  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
446  else
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;
450  }
451 
452  m_isInitialized = true;
453  bool ok = Traits::inplace_decomposition(m_matrix);
454  m_info = ok ? Success : NumericalIssue;
455 
456  return *this;
457 }
458 
464 template<typename MatrixType_, int UpLo_>
465 template<typename VectorType>
466 LLT<MatrixType_,UpLo_> & LLT<MatrixType_,UpLo_>::rankUpdate(const VectorType& v, const RealScalar& sigma)
467 {
469  eigen_assert(v.size()==m_matrix.cols());
471  if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
473  else
474  m_info = Success;
475 
476  return *this;
477 }
478 
479 #ifndef EIGEN_PARSED_BY_DOXYGEN
480 template<typename MatrixType_,int UpLo_>
481 template<typename RhsType, typename DstType>
482 void LLT<MatrixType_,UpLo_>::_solve_impl(const RhsType &rhs, DstType &dst) const
483 {
484  _solve_impl_transposed<true>(rhs, dst);
485 }
486 
487 template<typename MatrixType_,int UpLo_>
488 template<bool Conjugate, typename RhsType, typename DstType>
489 void LLT<MatrixType_,UpLo_>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
490 {
491  dst = rhs;
492 
493  matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
494  matrixU().template conjugateIf<!Conjugate>().solveInPlace(dst);
495 }
496 #endif
497 
511 template<typename MatrixType, int UpLo_>
512 template<typename Derived>
514 {
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);
519 }
520 
524 template<typename MatrixType, int UpLo_>
526 {
527  eigen_assert(m_isInitialized && "LLT is not initialized.");
528  return matrixL() * matrixL().adjoint().toDenseMatrix();
529 }
530 
535 template<typename Derived>
538 {
539  return LLT<PlainObject>(derived());
540 }
541 
546 template<typename MatrixType, unsigned int UpLo>
549 {
551 }
552 
553 } // end namespace Eigen
554 
555 #endif // EIGEN_LLT_H
Matrix3f m
const SqrtReturnType sqrt() const
Array< int, Dynamic, 1 > v
Array< int, 3, 1 > b
int n
ColXpr col(Index i)
This is the const version of col().
RealReturnType real() const
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
Definition: Macros.h:1149
#define EIGEN_NOEXCEPT
Definition: Macros.h:1260
#define EIGEN_CONSTEXPR
Definition: Macros.h:747
#define eigen_assert(x)
Definition: Macros.h:902
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:81
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:36
Matrix< float, 1, Dynamic > MatrixType
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:107
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:37
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:164
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
LLT & compute(const EigenBase< InputType > &matrix)
RealScalar m_l1_norm
Definition: LLT.h:229
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LLT.h:193
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: LLT.h:207
MatrixType m_matrix
Definition: LLT.h:228
LLT(Index size)
Default Constructor with memory preallocation.
Definition: LLT.h:103
LLT & rankUpdate(const VectorType &vec, const RealScalar &sigma=1)
@ MaxColsAtCompileTime
Definition: LLT.h:78
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: LLT.h:206
RealScalar rcond() const
Definition: LLT.h:168
bool m_isInitialized
Definition: LLT.h:230
MatrixType_ MatrixType
Definition: LLT.h:72
@ UpLo
Definition: LLT.h:84
@ PacketSize
Definition: LLT.h:82
@ AlignmentMask
Definition: LLT.h:83
LLT(EigenBase< InputType > &matrix)
Constructs a LLT factorization from a given matrix.
Definition: LLT.h:122
const LLT & adjoint() const EIGEN_NOEXCEPT
Definition: LLT.h:204
internal::LLT_Traits< MatrixType, UpLo > Traits
Definition: LLT.h:87
void solveInPlace(const MatrixBase< Derived > &bAndX) const
Definition: LLT.h:513
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
const MatrixType & matrixLLT() const
Definition: LLT.h:179
Traits::MatrixU matrixU() const
Definition: LLT.h:130
LLT(const EigenBase< InputType > &matrix)
Definition: LLT.h:107
LLT()
Default Constructor.
Definition: LLT.h:95
MatrixType reconstructedMatrix() const
Definition: LLT.h:525
ComputationInfo m_info
Definition: LLT.h:231
SolverBase< LLT > Base
Definition: LLT.h:73
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
internal::traits< LLT< MatrixType_, UpLo_ > >::Scalar Scalar
Definition: SolverBase.h:75
LLT< MatrixType_, UpLo_ > & derived()
Definition: EigenBase.h:48
ComputationInfo
Definition: Constants.h:444
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
static Index llt_rank_update_lower(MatrixType &mat, const VectorType &vec, const typename MatrixType::RealScalar &sigma)
Definition: LLT.h:239
const Scalar & y
void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
Definition: Jacobi.h:455
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
bool is_same_dense(const T1 &mat1, const T2 &mat2, std::enable_if_t< possibly_same_dense< T1, T2 >::value > *=0)
Definition: XprHelper.h:745
bool abs2(bool x)
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Derived & derived()
Definition: EigenBase.h:48
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
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:231
std::ptrdiff_t j