Eigen  3.2.10
 
Loading...
Searching...
No Matches
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
13namespace Eigen {
14
15namespace internal{
16template<typename MatrixType, int UpLo> struct LLT_Traits;
17}
18
46 /* HEY THIS DOX IS DISABLED BECAUSE THERE's A BUG EITHER HERE OR IN LDLT ABOUT THAT (OR BOTH)
47 * Note that during the decomposition, only the upper triangular part of A is considered. Therefore,
48 * the strict lower part does not have to store correct values.
49 */
50template<typename _MatrixType, int _UpLo> class LLT
51{
52 public:
53 typedef _MatrixType MatrixType;
54 enum {
55 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57 Options = MatrixType::Options,
58 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
59 };
60 typedef typename MatrixType::Scalar Scalar;
61 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
62 typedef typename MatrixType::Index Index;
63
64 enum {
65 PacketSize = internal::packet_traits<Scalar>::size,
66 AlignmentMask = int(PacketSize)-1,
67 UpLo = _UpLo
68 };
69
70 typedef internal::LLT_Traits<MatrixType,UpLo> Traits;
71
78 LLT() : m_matrix(), m_isInitialized(false) {}
79
86 LLT(Index size) : m_matrix(size, size),
87 m_isInitialized(false) {}
88
89 LLT(const MatrixType& matrix)
90 : m_matrix(matrix.rows(), matrix.cols()),
91 m_isInitialized(false)
92 {
93 compute(matrix);
94 }
95
97 inline typename Traits::MatrixU matrixU() const
98 {
99 eigen_assert(m_isInitialized && "LLT is not initialized.");
100 return Traits::getU(m_matrix);
101 }
102
104 inline typename Traits::MatrixL matrixL() const
105 {
106 eigen_assert(m_isInitialized && "LLT is not initialized.");
107 return Traits::getL(m_matrix);
108 }
109
120 template<typename Rhs>
121 inline const internal::solve_retval<LLT, Rhs>
122 solve(const MatrixBase<Rhs>& b) const
123 {
124 eigen_assert(m_isInitialized && "LLT is not initialized.");
125 eigen_assert(m_matrix.rows()==b.rows()
126 && "LLT::solve(): invalid number of rows of the right hand side matrix b");
127 return internal::solve_retval<LLT, Rhs>(*this, b.derived());
128 }
129
130 #ifdef EIGEN2_SUPPORT
131 template<typename OtherDerived, typename ResultType>
132 bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
133 {
134 *result = this->solve(b);
135 return true;
136 }
137
138 bool isPositiveDefinite() const { return true; }
139 #endif
140
141 template<typename Derived>
142 void solveInPlace(MatrixBase<Derived> &bAndX) const;
143
144 LLT& compute(const MatrixType& matrix);
145
150 inline const MatrixType& matrixLLT() const
151 {
152 eigen_assert(m_isInitialized && "LLT is not initialized.");
153 return m_matrix;
154 }
155
156 MatrixType reconstructedMatrix() const;
157
158
165 {
166 eigen_assert(m_isInitialized && "LLT is not initialized.");
167 return m_info;
168 }
169
170 inline Index rows() const { return m_matrix.rows(); }
171 inline Index cols() const { return m_matrix.cols(); }
172
173 template<typename VectorType>
174 LLT rankUpdate(const VectorType& vec, const RealScalar& sigma = 1);
175
176 protected:
177
178 static void check_template_parameters()
179 {
180 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
181 }
182
187 MatrixType m_matrix;
188 bool m_isInitialized;
189 ComputationInfo m_info;
190};
191
192namespace internal {
193
194template<typename Scalar, int UpLo> struct llt_inplace;
195
196template<typename MatrixType, typename VectorType>
197static typename MatrixType::Index llt_rank_update_lower(MatrixType& mat, const VectorType& vec, const typename MatrixType::RealScalar& sigma)
198{
199 using std::sqrt;
200 typedef typename MatrixType::Scalar Scalar;
201 typedef typename MatrixType::RealScalar RealScalar;
202 typedef typename MatrixType::Index Index;
203 typedef typename MatrixType::ColXpr ColXpr;
204 typedef typename internal::remove_all<ColXpr>::type ColXprCleaned;
205 typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
206 typedef Matrix<Scalar,Dynamic,1> TempVectorType;
207 typedef typename TempVectorType::SegmentReturnType TempVecSegment;
208
209 Index n = mat.cols();
210 eigen_assert(mat.rows()==n && vec.size()==n);
211
212 TempVectorType temp;
213
214 if(sigma>0)
215 {
216 // This version is based on Givens rotations.
217 // It is faster than the other one below, but only works for updates,
218 // i.e., for sigma > 0
219 temp = sqrt(sigma) * vec;
220
221 for(Index i=0; i<n; ++i)
222 {
223 JacobiRotation<Scalar> g;
224 g.makeGivens(mat(i,i), -temp(i), &mat(i,i));
225
226 Index rs = n-i-1;
227 if(rs>0)
228 {
229 ColXprSegment x(mat.col(i).tail(rs));
230 TempVecSegment y(temp.tail(rs));
231 apply_rotation_in_the_plane(x, y, g);
232 }
233 }
234 }
235 else
236 {
237 temp = vec;
238 RealScalar beta = 1;
239 for(Index j=0; j<n; ++j)
240 {
241 RealScalar Ljj = numext::real(mat.coeff(j,j));
242 RealScalar dj = numext::abs2(Ljj);
243 Scalar wj = temp.coeff(j);
244 RealScalar swj2 = sigma*numext::abs2(wj);
245 RealScalar gamma = dj*beta + swj2;
246
247 RealScalar x = dj + swj2/beta;
248 if (x<=RealScalar(0))
249 return j;
250 RealScalar nLjj = sqrt(x);
251 mat.coeffRef(j,j) = nLjj;
252 beta += swj2/dj;
253
254 // Update the terms of L
255 Index rs = n-j-1;
256 if(rs)
257 {
258 temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
259 if(gamma != 0)
260 mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
261 }
262 }
263 }
264 return -1;
265}
266
267template<typename Scalar> struct llt_inplace<Scalar, Lower>
268{
269 typedef typename NumTraits<Scalar>::Real RealScalar;
270 template<typename MatrixType>
271 static typename MatrixType::Index unblocked(MatrixType& mat)
272 {
273 using std::sqrt;
274 typedef typename MatrixType::Index Index;
275
276 eigen_assert(mat.rows()==mat.cols());
277 const Index size = mat.rows();
278 for(Index k = 0; k < size; ++k)
279 {
280 Index rs = size-k-1; // remaining size
281
282 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
283 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
284 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
285
286 RealScalar x = numext::real(mat.coeff(k,k));
287 if (k>0) x -= A10.squaredNorm();
288 if (x<=RealScalar(0))
289 return k;
290 mat.coeffRef(k,k) = x = sqrt(x);
291 if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
292 if (rs>0) A21 /= x;
293 }
294 return -1;
295 }
296
297 template<typename MatrixType>
298 static typename MatrixType::Index blocked(MatrixType& m)
299 {
300 typedef typename MatrixType::Index Index;
301 eigen_assert(m.rows()==m.cols());
302 Index size = m.rows();
303 if(size<32)
304 return unblocked(m);
305
306 Index blockSize = size/8;
307 blockSize = (blockSize/16)*16;
308 blockSize = (std::min)((std::max)(blockSize,Index(8)), Index(128));
309
310 for (Index k=0; k<size; k+=blockSize)
311 {
312 // partition the matrix:
313 // A00 | - | -
314 // lu = A10 | A11 | -
315 // A20 | A21 | A22
316 Index bs = (std::min)(blockSize, size-k);
317 Index rs = size - k - bs;
318 Block<MatrixType,Dynamic,Dynamic> A11(m,k, k, bs,bs);
319 Block<MatrixType,Dynamic,Dynamic> A21(m,k+bs,k, rs,bs);
320 Block<MatrixType,Dynamic,Dynamic> A22(m,k+bs,k+bs,rs,rs);
321
322 Index ret;
323 if((ret=unblocked(A11))>=0) return k+ret;
324 if(rs>0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
325 if(rs>0) A22.template selfadjointView<Lower>().rankUpdate(A21,-1); // bottleneck
326 }
327 return -1;
328 }
329
330 template<typename MatrixType, typename VectorType>
331 static typename MatrixType::Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
332 {
333 return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
334 }
335};
336
337template<typename Scalar> struct llt_inplace<Scalar, Upper>
338{
339 typedef typename NumTraits<Scalar>::Real RealScalar;
340
341 template<typename MatrixType>
342 static EIGEN_STRONG_INLINE typename MatrixType::Index unblocked(MatrixType& mat)
343 {
344 Transpose<MatrixType> matt(mat);
345 return llt_inplace<Scalar, Lower>::unblocked(matt);
346 }
347 template<typename MatrixType>
348 static EIGEN_STRONG_INLINE typename MatrixType::Index blocked(MatrixType& mat)
349 {
350 Transpose<MatrixType> matt(mat);
351 return llt_inplace<Scalar, Lower>::blocked(matt);
352 }
353 template<typename MatrixType, typename VectorType>
354 static typename MatrixType::Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
355 {
356 Transpose<MatrixType> matt(mat);
357 return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
358 }
359};
360
361template<typename MatrixType> struct LLT_Traits<MatrixType,Lower>
362{
363 typedef const TriangularView<const MatrixType, Lower> MatrixL;
365 static inline MatrixL getL(const MatrixType& m) { return m; }
366 static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
367 static bool inplace_decomposition(MatrixType& m)
368 { return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m)==-1; }
369};
370
371template<typename MatrixType> struct LLT_Traits<MatrixType,Upper>
372{
374 typedef const TriangularView<const MatrixType, Upper> MatrixU;
375 static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
376 static inline MatrixU getU(const MatrixType& m) { return m; }
377 static bool inplace_decomposition(MatrixType& m)
378 { return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m)==-1; }
379};
380
381} // end namespace internal
382
390template<typename MatrixType, int _UpLo>
392{
393 check_template_parameters();
394
395 eigen_assert(a.rows()==a.cols());
396 const Index size = a.rows();
397 m_matrix.resize(size, size);
398 m_matrix = a;
399
400 m_isInitialized = true;
401 bool ok = Traits::inplace_decomposition(m_matrix);
402 m_info = ok ? Success : NumericalIssue;
403
404 return *this;
405}
406
412template<typename _MatrixType, int _UpLo>
413template<typename VectorType>
414LLT<_MatrixType,_UpLo> LLT<_MatrixType,_UpLo>::rankUpdate(const VectorType& v, const RealScalar& sigma)
415{
416 EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
417 eigen_assert(v.size()==m_matrix.cols());
418 eigen_assert(m_isInitialized);
419 if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
420 m_info = NumericalIssue;
421 else
422 m_info = Success;
423
424 return *this;
425}
426
427namespace internal {
428template<typename _MatrixType, int UpLo, typename Rhs>
429struct solve_retval<LLT<_MatrixType, UpLo>, Rhs>
430 : solve_retval_base<LLT<_MatrixType, UpLo>, Rhs>
431{
432 typedef LLT<_MatrixType,UpLo> LLTType;
433 EIGEN_MAKE_SOLVE_HELPERS(LLTType,Rhs)
434
435 template<typename Dest> void evalTo(Dest& dst) const
436 {
437 dst = rhs();
438 dec().solveInPlace(dst);
439 }
440};
441}
442
456template<typename MatrixType, int _UpLo>
457template<typename Derived>
458void LLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
459{
460 eigen_assert(m_isInitialized && "LLT is not initialized.");
461 eigen_assert(m_matrix.rows()==bAndX.rows());
462 matrixL().solveInPlace(bAndX);
463 matrixU().solveInPlace(bAndX);
464}
465
469template<typename MatrixType, int _UpLo>
471{
472 eigen_assert(m_isInitialized && "LLT is not initialized.");
473 return matrixL() * matrixL().adjoint().toDenseMatrix();
474}
475
479template<typename Derived>
482{
483 return LLT<PlainObject>(derived());
484}
485
489template<typename MatrixType, unsigned int UpLo>
492{
493 return LLT<PlainObject,UpLo>(m_matrix);
494}
495
496} // end namespace Eigen
497
498#endif // EIGEN_LLT_H
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition LLT.h:51
Traits::MatrixL matrixL() const
Definition LLT.h:104
MatrixType reconstructedMatrix() const
Definition LLT.h:470
LLT()
Default Constructor.
Definition LLT.h:78
const MatrixType & matrixLLT() const
Definition LLT.h:150
Traits::MatrixU matrixU() const
Definition LLT.h:97
ComputationInfo info() const
Reports whether previous computation was successful.
Definition LLT.h:164
LLT & compute(const MatrixType &matrix)
Definition LLT.h:391
const internal::solve_retval< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition LLT.h:122
LLT(Index size)
Default Constructor with memory preallocation.
Definition LLT.h:86
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
const LLT< PlainObject > llt() const
Definition LLT.h:481
const LLT< PlainObject, UpLo > llt() const
Definition LLT.h:491
Base class for triangular part in a matrix.
Definition TriangularMatrix.h:160
ComputationInfo
Definition Constants.h:374
@ NumericalIssue
Definition Constants.h:378
@ Success
Definition Constants.h:376
@ Upper
Definition Constants.h:169
@ Lower
Definition Constants.h:167
Definition LDLT.h:18