$darkmode
Eigen  5.0.1-dev
LLT.h
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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
20 template <typename MatrixType_, int UpLo_>
21 struct traits<LLT<MatrixType_, UpLo_> > : traits<MatrixType_> {
22  typedef MatrixXpr XprKind;
23  typedef SolverStorage StorageKind;
24  typedef int StorageIndex;
25  enum { Flags = 0 };
26 };
27 
28 template <typename MatrixType, int UpLo>
29 struct LLT_Traits;
30 } // namespace internal
31 
69 template <typename MatrixType_, int UpLo_>
70 class LLT : public SolverBase<LLT<MatrixType_, UpLo_> > {
71  public:
72  typedef MatrixType_ MatrixType;
73  typedef SolverBase<LLT> Base;
74  friend class SolverBase<LLT>;
75 
76  EIGEN_GENERIC_PUBLIC_INTERFACE(LLT)
77  enum { MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
78 
79  enum { PacketSize = internal::packet_traits<Scalar>::size, AlignmentMask = int(PacketSize) - 1, UpLo = UpLo_ };
80 
81  typedef internal::LLT_Traits<MatrixType, UpLo> Traits;
82 
89  LLT() : m_matrix(), m_isInitialized(false) {}
90 
97  explicit LLT(Index size) : m_matrix(size, size), m_isInitialized(false) {}
98 
99  template <typename InputType>
100  explicit LLT(const EigenBase<InputType>& matrix) : m_matrix(matrix.rows(), matrix.cols()), m_isInitialized(false) {
101  compute(matrix.derived());
102  }
103 
111  template <typename InputType>
112  explicit LLT(EigenBase<InputType>& matrix) : m_matrix(matrix.derived()), m_isInitialized(false) {
113  compute(matrix.derived());
114  }
115 
117  inline typename Traits::MatrixU matrixU() const {
118  eigen_assert(m_isInitialized && "LLT is not initialized.");
119  return Traits::getU(m_matrix);
120  }
121 
123  inline typename Traits::MatrixL matrixL() const {
124  eigen_assert(m_isInitialized && "LLT is not initialized.");
125  return Traits::getL(m_matrix);
126  }
127 
128 #ifdef EIGEN_PARSED_BY_DOXYGEN
129 
139  template <typename Rhs>
140  inline const Solve<LLT, Rhs> solve(const MatrixBase<Rhs>& b) const;
141 #endif
142 
143  template <typename Derived>
144  void solveInPlace(const MatrixBase<Derived>& bAndX) const;
145 
146  template <typename InputType>
147  LLT& compute(const EigenBase<InputType>& matrix);
148 
152  RealScalar rcond() const {
153  eigen_assert(m_isInitialized && "LLT is not initialized.");
154  eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
155  return internal::rcond_estimate_helper(m_l1_norm, *this);
156  }
157 
162  inline const MatrixType& matrixLLT() const {
163  eigen_assert(m_isInitialized && "LLT is not initialized.");
164  return m_matrix;
165  }
166 
167  MatrixType reconstructedMatrix() const;
168 
175  eigen_assert(m_isInitialized && "LLT is not initialized.");
176  return m_info;
177  }
178 
185  const LLT& adjoint() const noexcept { return *this; }
186 
187  constexpr Index rows() const noexcept { return m_matrix.rows(); }
188  constexpr Index cols() const noexcept { return m_matrix.cols(); }
189 
190  template <typename VectorType>
191  LLT& rankUpdate(const VectorType& vec, const RealScalar& sigma = 1);
192 
193 #ifndef EIGEN_PARSED_BY_DOXYGEN
194  template <typename RhsType, typename DstType>
195  void _solve_impl(const RhsType& rhs, DstType& dst) const;
196 
197  template <bool Conjugate, typename RhsType, typename DstType>
198  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
199 #endif
200 
201  protected:
202  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
203 
204 
208  MatrixType m_matrix;
209  RealScalar m_l1_norm;
210  bool m_isInitialized;
211  ComputationInfo m_info;
212 };
213 
214 namespace internal {
215 
216 template <typename Scalar, int UpLo>
217 struct llt_inplace;
218 
219 template <typename MatrixType, typename VectorType>
220 static Index llt_rank_update_lower(MatrixType& mat, const VectorType& vec,
221  const typename MatrixType::RealScalar& sigma) {
222  using std::sqrt;
223  typedef typename MatrixType::Scalar Scalar;
224  typedef typename MatrixType::RealScalar RealScalar;
225  typedef typename MatrixType::ColXpr ColXpr;
226  typedef internal::remove_all_t<ColXpr> ColXprCleaned;
227  typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
228  typedef Matrix<Scalar, Dynamic, 1> TempVectorType;
229  typedef typename TempVectorType::SegmentReturnType TempVecSegment;
230 
231  Index n = mat.cols();
232  eigen_assert(mat.rows() == n && vec.size() == n);
233 
234  TempVectorType temp;
235 
236  if (sigma > 0) {
237  // This version is based on Givens rotations.
238  // It is faster than the other one below, but only works for updates,
239  // i.e., for sigma > 0
240  temp = sqrt(sigma) * vec;
241 
242  for (Index i = 0; i < n; ++i) {
243  JacobiRotation<Scalar> g;
244  g.makeGivens(mat(i, i), -temp(i), &mat(i, i));
245 
246  Index rs = n - i - 1;
247  if (rs > 0) {
248  ColXprSegment x(mat.col(i).tail(rs));
249  TempVecSegment y(temp.tail(rs));
250  apply_rotation_in_the_plane(x, y, g);
251  }
252  }
253  } else {
254  temp = vec;
255  RealScalar beta = 1;
256  for (Index j = 0; j < n; ++j) {
257  RealScalar Ljj = numext::real(mat.coeff(j, j));
258  RealScalar dj = numext::abs2(Ljj);
259  Scalar wj = temp.coeff(j);
260  RealScalar swj2 = sigma * numext::abs2(wj);
261  RealScalar gamma = dj * beta + swj2;
262 
263  RealScalar x = dj + swj2 / beta;
264  if (x <= RealScalar(0)) return j;
265  RealScalar nLjj = sqrt(x);
266  mat.coeffRef(j, j) = nLjj;
267  beta += swj2 / dj;
268 
269  // Update the terms of L
270  Index rs = n - j - 1;
271  if (rs) {
272  temp.tail(rs) -= (wj / Ljj) * mat.col(j).tail(rs);
273  if (!numext::is_exactly_zero(gamma))
274  mat.col(j).tail(rs) =
275  (nLjj / Ljj) * mat.col(j).tail(rs) + (nLjj * sigma * numext::conj(wj) / gamma) * temp.tail(rs);
276  }
277  }
278  }
279  return -1;
280 }
281 
282 template <typename Scalar>
283 struct llt_inplace<Scalar, Lower> {
284  typedef typename NumTraits<Scalar>::Real RealScalar;
285  template <typename MatrixType>
286  static Index unblocked(MatrixType& mat) {
287  using std::sqrt;
288 
289  eigen_assert(mat.rows() == mat.cols());
290  const Index size = mat.rows();
291  for (Index k = 0; k < size; ++k) {
292  Index rs = size - k - 1; // remaining size
293 
294  Block<MatrixType, Dynamic, 1> A21(mat, k + 1, k, rs, 1);
295  Block<MatrixType, 1, Dynamic> A10(mat, k, 0, 1, k);
296  Block<MatrixType, Dynamic, Dynamic> A20(mat, k + 1, 0, rs, k);
297 
298  RealScalar x = numext::real(mat.coeff(k, k));
299  if (k > 0) x -= A10.squaredNorm();
300  if (x <= RealScalar(0)) return k;
301  mat.coeffRef(k, k) = x = sqrt(x);
302  if (k > 0 && rs > 0) A21.noalias() -= A20 * A10.adjoint();
303  if (rs > 0) A21 /= x;
304  }
305  return -1;
306  }
307 
308  template <typename MatrixType>
309  static Index blocked(MatrixType& m) {
310  eigen_assert(m.rows() == m.cols());
311  Index size = m.rows();
312  if (size < 32) return unblocked(m);
313 
314  Index blockSize = size / 8;
315  blockSize = (blockSize / 16) * 16;
316  blockSize = (std::min)((std::max)(blockSize, Index(8)), Index(128));
317 
318  for (Index k = 0; k < size; k += blockSize) {
319  // partition the matrix:
320  // A00 | - | -
321  // lu = A10 | A11 | -
322  // A20 | A21 | A22
323  Index bs = (std::min)(blockSize, size - k);
324  Index rs = size - k - bs;
325  Block<MatrixType, Dynamic, Dynamic> A11(m, k, k, bs, bs);
326  Block<MatrixType, Dynamic, Dynamic> A21(m, k + bs, k, rs, bs);
327  Block<MatrixType, Dynamic, Dynamic> A22(m, k + bs, k + bs, rs, rs);
328 
329  Index ret;
330  if ((ret = unblocked(A11)) >= 0) return k + ret;
331  if (rs > 0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
332  if (rs > 0)
333  A22.template selfadjointView<Lower>().rankUpdate(A21,
334  typename NumTraits<RealScalar>::Literal(-1)); // bottleneck
335  }
336  return -1;
337  }
338 
339  template <typename MatrixType, typename VectorType>
340  static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma) {
341  return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
342  }
343 };
344 
345 template <typename Scalar>
346 struct llt_inplace<Scalar, Upper> {
347  typedef typename NumTraits<Scalar>::Real RealScalar;
348 
349  template <typename MatrixType>
350  static EIGEN_STRONG_INLINE Index unblocked(MatrixType& mat) {
351  Transpose<MatrixType> matt(mat);
352  return llt_inplace<Scalar, Lower>::unblocked(matt);
353  }
354  template <typename MatrixType>
355  static EIGEN_STRONG_INLINE Index blocked(MatrixType& mat) {
356  Transpose<MatrixType> matt(mat);
357  return llt_inplace<Scalar, Lower>::blocked(matt);
358  }
359  template <typename MatrixType, typename VectorType>
360  static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma) {
361  Transpose<MatrixType> matt(mat);
362  return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
363  }
364 };
365 
366 template <typename MatrixType>
367 struct LLT_Traits<MatrixType, Lower> {
368  typedef const TriangularView<const MatrixType, Lower> MatrixL;
369  typedef const TriangularView<const typename MatrixType::AdjointReturnType, Upper> MatrixU;
370  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
371  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
372  static bool inplace_decomposition(MatrixType& m) {
373  return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m) == -1;
374  }
375 };
376 
377 template <typename MatrixType>
378 struct LLT_Traits<MatrixType, Upper> {
379  typedef const TriangularView<const typename MatrixType::AdjointReturnType, Lower> MatrixL;
380  typedef const TriangularView<const MatrixType, Upper> MatrixU;
381  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
382  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
383  static bool inplace_decomposition(MatrixType& m) {
384  return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m) == -1;
385  }
386 };
387 
388 } // end namespace internal
389 
397 template <typename MatrixType, int UpLo_>
398 template <typename InputType>
400  eigen_assert(a.rows() == a.cols());
401  const Index size = a.rows();
402  m_matrix.resize(size, size);
403  if (!internal::is_same_dense(m_matrix, a.derived())) m_matrix = a.derived();
404 
405  // Compute matrix L1 norm = max abs column sum.
406  m_l1_norm = RealScalar(0);
407  // TODO move this code to SelfAdjointView
408  for (Index col = 0; col < size; ++col) {
409  RealScalar abs_col_sum;
410  if (UpLo_ == Lower)
411  abs_col_sum =
412  m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
413  else
414  abs_col_sum =
415  m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
416  if (abs_col_sum > m_l1_norm) m_l1_norm = abs_col_sum;
417  }
418 
419  m_isInitialized = true;
420  bool ok = Traits::inplace_decomposition(m_matrix);
421  m_info = ok ? Success : NumericalIssue;
422 
423  return *this;
424 }
425 
431 template <typename MatrixType_, int UpLo_>
432 template <typename VectorType>
433 LLT<MatrixType_, UpLo_>& LLT<MatrixType_, UpLo_>::rankUpdate(const VectorType& v, const RealScalar& sigma) {
434  EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
435  eigen_assert(v.size() == m_matrix.cols());
436  eigen_assert(m_isInitialized);
437  if (internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix, v, sigma) >= 0)
438  m_info = NumericalIssue;
439  else
440  m_info = Success;
441 
442  return *this;
443 }
444 
445 #ifndef EIGEN_PARSED_BY_DOXYGEN
446 template <typename MatrixType_, int UpLo_>
447 template <typename RhsType, typename DstType>
448 void LLT<MatrixType_, UpLo_>::_solve_impl(const RhsType& rhs, DstType& dst) const {
449  _solve_impl_transposed<true>(rhs, dst);
450 }
451 
452 template <typename MatrixType_, int UpLo_>
453 template <bool Conjugate, typename RhsType, typename DstType>
454 void LLT<MatrixType_, UpLo_>::_solve_impl_transposed(const RhsType& rhs, DstType& dst) const {
455  dst = rhs;
456 
457  matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
458  matrixU().template conjugateIf<!Conjugate>().solveInPlace(dst);
459 }
460 #endif
461 
475 template <typename MatrixType, int UpLo_>
476 template <typename Derived>
477 void LLT<MatrixType, UpLo_>::solveInPlace(const MatrixBase<Derived>& bAndX) const {
478  eigen_assert(m_isInitialized && "LLT is not initialized.");
479  eigen_assert(m_matrix.rows() == bAndX.rows());
480  matrixL().solveInPlace(bAndX);
481  matrixU().solveInPlace(bAndX);
482 }
483 
487 template <typename MatrixType, int UpLo_>
489  eigen_assert(m_isInitialized && "LLT is not initialized.");
490  return matrixL() * matrixL().adjoint().toDenseMatrix();
491 }
492 
497 template <typename Derived>
499  return LLT<PlainObject>(derived());
500 }
501 
506 template <typename MatrixType, unsigned int UpLo>
508  const {
509  return LLT<PlainObject, UpLo>(m_matrix);
510 }
511 
512 } // end namespace Eigen
513 
514 #endif // EIGEN_LLT_H
constexpr Index size() const noexcept
Definition: EigenBase.h:64
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
RealScalar rcond() const
Definition: LLT.h:152
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
MatrixType reconstructedMatrix() const
Definition: LLT.h:488
constexpr Index rows() const noexcept
Definition: EigenBase.h:59
Definition: EigenBase.h:33
constexpr Index cols() const noexcept
Definition: EigenBase.h:61
Definition: Constants.h:211
LLT(EigenBase< InputType > &matrix)
Constructs a LLT factorization from a given matrix.
Definition: LLT.h:112
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
Definition: Constants.h:442
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
LLT(Index size)
Default Constructor with memory preallocation.
Definition: LLT.h:97
const MatrixType & matrixLLT() const
Definition: LLT.h:162
const LLT & adjoint() const noexcept
Definition: LLT.h:185
Definition: Constants.h:213
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
LLT()
Default Constructor.
Definition: LLT.h:89
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LLT.h:174
Traits::MatrixU matrixU() const
Definition: LLT.h:117
const LLT< PlainObject > llt() const
Definition: LLT.h:498
Pseudo expression representing a solving operation.
Definition: Solve.h:62
ComputationInfo
Definition: Constants.h:438
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:72
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:507
Traits::MatrixL matrixL() const
Definition: LLT.h:123