• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 namespace Eigen {
14 
15 namespace internal{
16 template<typename MatrixType, int UpLo> struct LLT_Traits;
17 }
18 
19 /** \ingroup Cholesky_Module
20   *
21   * \class LLT
22   *
23   * \brief Standard Cholesky decomposition (LL^T) of a matrix and associated features
24   *
25   * \tparam _MatrixType the type of the matrix of which we are computing the LL^T Cholesky decomposition
26   * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
27   *             The other triangular part won't be read.
28   *
29   * This class performs a LL^T Cholesky decomposition of a symmetric, positive definite
30   * matrix A such that A = LL^* = U^*U, where L is lower triangular.
31   *
32   * While the Cholesky decomposition is particularly useful to solve selfadjoint problems like  D^*D x = b,
33   * for that purpose, we recommend the Cholesky decomposition without square root which is more stable
34   * and even faster. Nevertheless, this standard Cholesky decomposition remains useful in many other
35   * situations like generalised eigen problems with hermitian matrices.
36   *
37   * Remember that Cholesky decompositions are not rank-revealing. This LLT decomposition is only stable on positive definite matrices,
38   * use LDLT instead for the semidefinite case. Also, do not use a Cholesky decomposition to determine whether a system of equations
39   * has a solution.
40   *
41   * Example: \include LLT_example.cpp
42   * Output: \verbinclude LLT_example.out
43   *
44   * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
45   *
46   * \sa MatrixBase::llt(), SelfAdjointView::llt(), class LDLT
47   */
48  /* HEY THIS DOX IS DISABLED BECAUSE THERE's A BUG EITHER HERE OR IN LDLT ABOUT THAT (OR BOTH)
49   * Note that during the decomposition, only the upper triangular part of A is considered. Therefore,
50   * the strict lower part does not have to store correct values.
51   */
52 template<typename _MatrixType, int _UpLo> class LLT
53 {
54   public:
55     typedef _MatrixType MatrixType;
56     enum {
57       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
58       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
59       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
60     };
61     typedef typename MatrixType::Scalar Scalar;
62     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
63     typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
64     typedef typename MatrixType::StorageIndex StorageIndex;
65 
66     enum {
67       PacketSize = internal::packet_traits<Scalar>::size,
68       AlignmentMask = int(PacketSize)-1,
69       UpLo = _UpLo
70     };
71 
72     typedef internal::LLT_Traits<MatrixType,UpLo> Traits;
73 
74     /**
75       * \brief Default Constructor.
76       *
77       * The default constructor is useful in cases in which the user intends to
78       * perform decompositions via LLT::compute(const MatrixType&).
79       */
LLT()80     LLT() : m_matrix(), m_isInitialized(false) {}
81 
82     /** \brief Default Constructor with memory preallocation
83       *
84       * Like the default constructor but with preallocation of the internal data
85       * according to the specified problem \a size.
86       * \sa LLT()
87       */
LLT(Index size)88     explicit LLT(Index size) : m_matrix(size, size),
89                     m_isInitialized(false) {}
90 
91     template<typename InputType>
LLT(const EigenBase<InputType> & matrix)92     explicit LLT(const EigenBase<InputType>& matrix)
93       : m_matrix(matrix.rows(), matrix.cols()),
94         m_isInitialized(false)
95     {
96       compute(matrix.derived());
97     }
98 
99     /** \brief Constructs a LDLT factorization from a given matrix
100       *
101       * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when
102       * \c MatrixType is a Eigen::Ref.
103       *
104       * \sa LLT(const EigenBase&)
105       */
106     template<typename InputType>
LLT(EigenBase<InputType> & matrix)107     explicit LLT(EigenBase<InputType>& matrix)
108       : m_matrix(matrix.derived()),
109         m_isInitialized(false)
110     {
111       compute(matrix.derived());
112     }
113 
114     /** \returns a view of the upper triangular matrix U */
matrixU()115     inline typename Traits::MatrixU matrixU() const
116     {
117       eigen_assert(m_isInitialized && "LLT is not initialized.");
118       return Traits::getU(m_matrix);
119     }
120 
121     /** \returns a view of the lower triangular matrix L */
matrixL()122     inline typename Traits::MatrixL matrixL() const
123     {
124       eigen_assert(m_isInitialized && "LLT is not initialized.");
125       return Traits::getL(m_matrix);
126     }
127 
128     /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A.
129       *
130       * Since this LLT class assumes anyway that the matrix A is invertible, the solution
131       * theoretically exists and is unique regardless of b.
132       *
133       * Example: \include LLT_solve.cpp
134       * Output: \verbinclude LLT_solve.out
135       *
136       * \sa solveInPlace(), MatrixBase::llt(), SelfAdjointView::llt()
137       */
138     template<typename Rhs>
139     inline const Solve<LLT, Rhs>
solve(const MatrixBase<Rhs> & b)140     solve(const MatrixBase<Rhs>& b) const
141     {
142       eigen_assert(m_isInitialized && "LLT is not initialized.");
143       eigen_assert(m_matrix.rows()==b.rows()
144                 && "LLT::solve(): invalid number of rows of the right hand side matrix b");
145       return Solve<LLT, Rhs>(*this, b.derived());
146     }
147 
148     template<typename Derived>
149     void solveInPlace(MatrixBase<Derived> &bAndX) const;
150 
151     template<typename InputType>
152     LLT& compute(const EigenBase<InputType>& matrix);
153 
154     /** \returns an estimate of the reciprocal condition number of the matrix of
155       *  which \c *this is the Cholesky decomposition.
156       */
rcond()157     RealScalar rcond() const
158     {
159       eigen_assert(m_isInitialized && "LLT is not initialized.");
160       eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
161       return internal::rcond_estimate_helper(m_l1_norm, *this);
162     }
163 
164     /** \returns the LLT decomposition matrix
165       *
166       * TODO: document the storage layout
167       */
matrixLLT()168     inline const MatrixType& matrixLLT() const
169     {
170       eigen_assert(m_isInitialized && "LLT is not initialized.");
171       return m_matrix;
172     }
173 
174     MatrixType reconstructedMatrix() const;
175 
176 
177     /** \brief Reports whether previous computation was successful.
178       *
179       * \returns \c Success if computation was succesful,
180       *          \c NumericalIssue if the matrix.appears to be negative.
181       */
info()182     ComputationInfo info() const
183     {
184       eigen_assert(m_isInitialized && "LLT is not initialized.");
185       return m_info;
186     }
187 
188     /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.
189       *
190       * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:
191       * \code x = decomposition.adjoint().solve(b) \endcode
192       */
adjoint()193     const LLT& adjoint() const { return *this; };
194 
rows()195     inline Index rows() const { return m_matrix.rows(); }
cols()196     inline Index cols() const { return m_matrix.cols(); }
197 
198     template<typename VectorType>
199     LLT rankUpdate(const VectorType& vec, const RealScalar& sigma = 1);
200 
201     #ifndef EIGEN_PARSED_BY_DOXYGEN
202     template<typename RhsType, typename DstType>
203     EIGEN_DEVICE_FUNC
204     void _solve_impl(const RhsType &rhs, DstType &dst) const;
205     #endif
206 
207   protected:
208 
check_template_parameters()209     static void check_template_parameters()
210     {
211       EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
212     }
213 
214     /** \internal
215       * Used to compute and store L
216       * The strict upper part is not used and even not initialized.
217       */
218     MatrixType m_matrix;
219     RealScalar m_l1_norm;
220     bool m_isInitialized;
221     ComputationInfo m_info;
222 };
223 
224 namespace internal {
225 
226 template<typename Scalar, int UpLo> struct llt_inplace;
227 
228 template<typename MatrixType, typename VectorType>
llt_rank_update_lower(MatrixType & mat,const VectorType & vec,const typename MatrixType::RealScalar & sigma)229 static Index llt_rank_update_lower(MatrixType& mat, const VectorType& vec, const typename MatrixType::RealScalar& sigma)
230 {
231   using std::sqrt;
232   typedef typename MatrixType::Scalar Scalar;
233   typedef typename MatrixType::RealScalar RealScalar;
234   typedef typename MatrixType::ColXpr ColXpr;
235   typedef typename internal::remove_all<ColXpr>::type ColXprCleaned;
236   typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
237   typedef Matrix<Scalar,Dynamic,1> TempVectorType;
238   typedef typename TempVectorType::SegmentReturnType TempVecSegment;
239 
240   Index n = mat.cols();
241   eigen_assert(mat.rows()==n && vec.size()==n);
242 
243   TempVectorType temp;
244 
245   if(sigma>0)
246   {
247     // This version is based on Givens rotations.
248     // It is faster than the other one below, but only works for updates,
249     // i.e., for sigma > 0
250     temp = sqrt(sigma) * vec;
251 
252     for(Index i=0; i<n; ++i)
253     {
254       JacobiRotation<Scalar> g;
255       g.makeGivens(mat(i,i), -temp(i), &mat(i,i));
256 
257       Index rs = n-i-1;
258       if(rs>0)
259       {
260         ColXprSegment x(mat.col(i).tail(rs));
261         TempVecSegment y(temp.tail(rs));
262         apply_rotation_in_the_plane(x, y, g);
263       }
264     }
265   }
266   else
267   {
268     temp = vec;
269     RealScalar beta = 1;
270     for(Index j=0; j<n; ++j)
271     {
272       RealScalar Ljj = numext::real(mat.coeff(j,j));
273       RealScalar dj = numext::abs2(Ljj);
274       Scalar wj = temp.coeff(j);
275       RealScalar swj2 = sigma*numext::abs2(wj);
276       RealScalar gamma = dj*beta + swj2;
277 
278       RealScalar x = dj + swj2/beta;
279       if (x<=RealScalar(0))
280         return j;
281       RealScalar nLjj = sqrt(x);
282       mat.coeffRef(j,j) = nLjj;
283       beta += swj2/dj;
284 
285       // Update the terms of L
286       Index rs = n-j-1;
287       if(rs)
288       {
289         temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
290         if(gamma != 0)
291           mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
292       }
293     }
294   }
295   return -1;
296 }
297 
298 template<typename Scalar> struct llt_inplace<Scalar, Lower>
299 {
300   typedef typename NumTraits<Scalar>::Real RealScalar;
301   template<typename MatrixType>
302   static Index unblocked(MatrixType& mat)
303   {
304     using std::sqrt;
305 
306     eigen_assert(mat.rows()==mat.cols());
307     const Index size = mat.rows();
308     for(Index k = 0; k < size; ++k)
309     {
310       Index rs = size-k-1; // remaining size
311 
312       Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
313       Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
314       Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
315 
316       RealScalar x = numext::real(mat.coeff(k,k));
317       if (k>0) x -= A10.squaredNorm();
318       if (x<=RealScalar(0))
319         return k;
320       mat.coeffRef(k,k) = x = sqrt(x);
321       if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
322       if (rs>0) A21 /= x;
323     }
324     return -1;
325   }
326 
327   template<typename MatrixType>
328   static Index blocked(MatrixType& m)
329   {
330     eigen_assert(m.rows()==m.cols());
331     Index size = m.rows();
332     if(size<32)
333       return unblocked(m);
334 
335     Index blockSize = size/8;
336     blockSize = (blockSize/16)*16;
337     blockSize = (std::min)((std::max)(blockSize,Index(8)), Index(128));
338 
339     for (Index k=0; k<size; k+=blockSize)
340     {
341       // partition the matrix:
342       //       A00 |  -  |  -
343       // lu  = A10 | A11 |  -
344       //       A20 | A21 | A22
345       Index bs = (std::min)(blockSize, size-k);
346       Index rs = size - k - bs;
347       Block<MatrixType,Dynamic,Dynamic> A11(m,k,   k,   bs,bs);
348       Block<MatrixType,Dynamic,Dynamic> A21(m,k+bs,k,   rs,bs);
349       Block<MatrixType,Dynamic,Dynamic> A22(m,k+bs,k+bs,rs,rs);
350 
351       Index ret;
352       if((ret=unblocked(A11))>=0) return k+ret;
353       if(rs>0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
354       if(rs>0) A22.template selfadjointView<Lower>().rankUpdate(A21,typename NumTraits<RealScalar>::Literal(-1)); // bottleneck
355     }
356     return -1;
357   }
358 
359   template<typename MatrixType, typename VectorType>
360   static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
361   {
362     return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
363   }
364 };
365 
366 template<typename Scalar> struct llt_inplace<Scalar, Upper>
367 {
368   typedef typename NumTraits<Scalar>::Real RealScalar;
369 
370   template<typename MatrixType>
371   static EIGEN_STRONG_INLINE Index unblocked(MatrixType& mat)
372   {
373     Transpose<MatrixType> matt(mat);
374     return llt_inplace<Scalar, Lower>::unblocked(matt);
375   }
376   template<typename MatrixType>
377   static EIGEN_STRONG_INLINE Index blocked(MatrixType& mat)
378   {
379     Transpose<MatrixType> matt(mat);
380     return llt_inplace<Scalar, Lower>::blocked(matt);
381   }
382   template<typename MatrixType, typename VectorType>
383   static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
384   {
385     Transpose<MatrixType> matt(mat);
386     return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
387   }
388 };
389 
390 template<typename MatrixType> struct LLT_Traits<MatrixType,Lower>
391 {
392   typedef const TriangularView<const MatrixType, Lower> MatrixL;
393   typedef const TriangularView<const typename MatrixType::AdjointReturnType, Upper> MatrixU;
394   static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
395   static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
396   static bool inplace_decomposition(MatrixType& m)
397   { return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m)==-1; }
398 };
399 
400 template<typename MatrixType> struct LLT_Traits<MatrixType,Upper>
401 {
402   typedef const TriangularView<const typename MatrixType::AdjointReturnType, Lower> MatrixL;
403   typedef const TriangularView<const MatrixType, Upper> MatrixU;
404   static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
405   static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
406   static bool inplace_decomposition(MatrixType& m)
407   { return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m)==-1; }
408 };
409 
410 } // end namespace internal
411 
412 /** Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of \a matrix
413   *
414   * \returns a reference to *this
415   *
416   * Example: \include TutorialLinAlgComputeTwice.cpp
417   * Output: \verbinclude TutorialLinAlgComputeTwice.out
418   */
419 template<typename MatrixType, int _UpLo>
420 template<typename InputType>
421 LLT<MatrixType,_UpLo>& LLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
422 {
423   check_template_parameters();
424 
425   eigen_assert(a.rows()==a.cols());
426   const Index size = a.rows();
427   m_matrix.resize(size, size);
428   m_matrix = a.derived();
429 
430   // Compute matrix L1 norm = max abs column sum.
431   m_l1_norm = RealScalar(0);
432   // TODO move this code to SelfAdjointView
433   for (Index col = 0; col < size; ++col) {
434     RealScalar abs_col_sum;
435     if (_UpLo == Lower)
436       abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
437     else
438       abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
439     if (abs_col_sum > m_l1_norm)
440       m_l1_norm = abs_col_sum;
441   }
442 
443   m_isInitialized = true;
444   bool ok = Traits::inplace_decomposition(m_matrix);
445   m_info = ok ? Success : NumericalIssue;
446 
447   return *this;
448 }
449 
450 /** Performs a rank one update (or dowdate) of the current decomposition.
451   * If A = LL^* before the rank one update,
452   * then after it we have LL^* = A + sigma * v v^* where \a v must be a vector
453   * of same dimension.
454   */
455 template<typename _MatrixType, int _UpLo>
456 template<typename VectorType>
457 LLT<_MatrixType,_UpLo> LLT<_MatrixType,_UpLo>::rankUpdate(const VectorType& v, const RealScalar& sigma)
458 {
459   EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
460   eigen_assert(v.size()==m_matrix.cols());
461   eigen_assert(m_isInitialized);
462   if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
463     m_info = NumericalIssue;
464   else
465     m_info = Success;
466 
467   return *this;
468 }
469 
470 #ifndef EIGEN_PARSED_BY_DOXYGEN
471 template<typename _MatrixType,int _UpLo>
472 template<typename RhsType, typename DstType>
473 void LLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
474 {
475   dst = rhs;
476   solveInPlace(dst);
477 }
478 #endif
479 
480 /** \internal use x = llt_object.solve(x);
481   *
482   * This is the \em in-place version of solve().
483   *
484   * \param bAndX represents both the right-hand side matrix b and result x.
485   *
486   * This version avoids a copy when the right hand side matrix b is not needed anymore.
487   *
488   * \sa LLT::solve(), MatrixBase::llt()
489   */
490 template<typename MatrixType, int _UpLo>
491 template<typename Derived>
492 void LLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
493 {
494   eigen_assert(m_isInitialized && "LLT is not initialized.");
495   eigen_assert(m_matrix.rows()==bAndX.rows());
496   matrixL().solveInPlace(bAndX);
497   matrixU().solveInPlace(bAndX);
498 }
499 
500 /** \returns the matrix represented by the decomposition,
501  * i.e., it returns the product: L L^*.
502  * This function is provided for debug purpose. */
503 template<typename MatrixType, int _UpLo>
504 MatrixType LLT<MatrixType,_UpLo>::reconstructedMatrix() const
505 {
506   eigen_assert(m_isInitialized && "LLT is not initialized.");
507   return matrixL() * matrixL().adjoint().toDenseMatrix();
508 }
509 
510 /** \cholesky_module
511   * \returns the LLT decomposition of \c *this
512   * \sa SelfAdjointView::llt()
513   */
514 template<typename Derived>
515 inline const LLT<typename MatrixBase<Derived>::PlainObject>
516 MatrixBase<Derived>::llt() const
517 {
518   return LLT<PlainObject>(derived());
519 }
520 
521 /** \cholesky_module
522   * \returns the LLT decomposition of \c *this
523   * \sa SelfAdjointView::llt()
524   */
525 template<typename MatrixType, unsigned int UpLo>
526 inline const LLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
527 SelfAdjointView<MatrixType, UpLo>::llt() const
528 {
529   return LLT<PlainObject,UpLo>(m_matrix);
530 }
531 
532 } // end namespace Eigen
533 
534 #endif // EIGEN_LLT_H
535