source: pacpussensors/trunk/Vislab/lib3dv/eigen/Eigen/src/Cholesky/LDLT.h@ 136

Last change on this file since 136 was 136, checked in by ldecherf, 7 years ago

Doc

File size: 21.2 KB
Line 
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
6// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7// Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
8//
9// This Source Code Form is subject to the terms of the Mozilla
10// Public License v. 2.0. If a copy of the MPL was not distributed
11// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12
13#ifndef EIGEN_LDLT_H
14#define EIGEN_LDLT_H
15
16namespace Eigen {
17
18namespace internal {
19 template<typename MatrixType, int UpLo> struct LDLT_Traits;
20
21 // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
22 enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
23}
24
25/** \ingroup Cholesky_Module
26 *
27 * \class LDLT
28 *
29 * \brief Robust Cholesky decomposition of a matrix with pivoting
30 *
31 * \param MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
32 * \param UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
33 * The other triangular part won't be read.
34 *
35 * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
36 * matrix \f$ A \f$ such that \f$ A = P^TLDL^*P \f$, where P is a permutation matrix, L
37 * is lower triangular with a unit diagonal and D is a diagonal matrix.
38 *
39 * The decomposition uses pivoting to ensure stability, so that L will have
40 * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
41 * on D also stabilizes the computation.
42 *
43 * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
44 * decomposition to determine whether a system of equations has a solution.
45 *
46 * \sa MatrixBase::ldlt(), class LLT
47 */
48template<typename _MatrixType, int _UpLo> class LDLT
49{
50 public:
51 typedef _MatrixType MatrixType;
52 enum {
53 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
54 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
55 Options = MatrixType::Options & ~RowMajorBit, // these are the options for the TmpMatrixType, we need a ColMajor matrix here!
56 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
57 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
58 UpLo = _UpLo
59 };
60 typedef typename MatrixType::Scalar Scalar;
61 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
62 typedef typename MatrixType::Index Index;
63 typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> TmpMatrixType;
64
65 typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
66 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
67
68 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
69
70 /** \brief Default Constructor.
71 *
72 * The default constructor is useful in cases in which the user intends to
73 * perform decompositions via LDLT::compute(const MatrixType&).
74 */
75 LDLT()
76 : m_matrix(),
77 m_transpositions(),
78 m_sign(internal::ZeroSign),
79 m_isInitialized(false)
80 {}
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 LDLT()
87 */
88 LDLT(Index size)
89 : m_matrix(size, size),
90 m_transpositions(size),
91 m_temporary(size),
92 m_sign(internal::ZeroSign),
93 m_isInitialized(false)
94 {}
95
96 /** \brief Constructor with decomposition
97 *
98 * This calculates the decomposition for the input \a matrix.
99 * \sa LDLT(Index size)
100 */
101 LDLT(const MatrixType& matrix)
102 : m_matrix(matrix.rows(), matrix.cols()),
103 m_transpositions(matrix.rows()),
104 m_temporary(matrix.rows()),
105 m_sign(internal::ZeroSign),
106 m_isInitialized(false)
107 {
108 compute(matrix);
109 }
110
111 /** Clear any existing decomposition
112 * \sa rankUpdate(w,sigma)
113 */
114 void setZero()
115 {
116 m_isInitialized = false;
117 }
118
119 /** \returns a view of the upper triangular matrix U */
120 inline typename Traits::MatrixU matrixU() const
121 {
122 eigen_assert(m_isInitialized && "LDLT is not initialized.");
123 return Traits::getU(m_matrix);
124 }
125
126 /** \returns a view of the lower triangular matrix L */
127 inline typename Traits::MatrixL matrixL() const
128 {
129 eigen_assert(m_isInitialized && "LDLT is not initialized.");
130 return Traits::getL(m_matrix);
131 }
132
133 /** \returns the permutation matrix P as a transposition sequence.
134 */
135 inline const TranspositionType& transpositionsP() const
136 {
137 eigen_assert(m_isInitialized && "LDLT is not initialized.");
138 return m_transpositions;
139 }
140
141 /** \returns the coefficients of the diagonal matrix D */
142 inline Diagonal<const MatrixType> vectorD() const
143 {
144 eigen_assert(m_isInitialized && "LDLT is not initialized.");
145 return m_matrix.diagonal();
146 }
147
148 /** \returns true if the matrix is positive (semidefinite) */
149 inline bool isPositive() const
150 {
151 eigen_assert(m_isInitialized && "LDLT is not initialized.");
152 return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
153 }
154
155 #ifdef EIGEN2_SUPPORT
156 inline bool isPositiveDefinite() const
157 {
158 return isPositive();
159 }
160 #endif
161
162 /** \returns true if the matrix is negative (semidefinite) */
163 inline bool isNegative(void) const
164 {
165 eigen_assert(m_isInitialized && "LDLT is not initialized.");
166 return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
167 }
168
169 /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
170 *
171 * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
172 *
173 * \note_about_checking_solutions
174 *
175 * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
176 * by solving the systems \f$ P^T y_1 = b \f$, \f$ L y_2 = y_1 \f$, \f$ D y_3 = y_2 \f$,
177 * \f$ L^* y_4 = y_3 \f$ and \f$ P x = y_4 \f$ in succession. If the matrix \f$ A \f$ is singular, then
178 * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
179 * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
180 * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
181 *
182 * \sa MatrixBase::ldlt()
183 */
184 template<typename Rhs>
185 inline const internal::solve_retval<LDLT, Rhs>
186 solve(const MatrixBase<Rhs>& b) const
187 {
188 eigen_assert(m_isInitialized && "LDLT is not initialized.");
189 eigen_assert(m_matrix.rows()==b.rows()
190 && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
191 return internal::solve_retval<LDLT, Rhs>(*this, b.derived());
192 }
193
194 #ifdef EIGEN2_SUPPORT
195 template<typename OtherDerived, typename ResultType>
196 bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
197 {
198 *result = this->solve(b);
199 return true;
200 }
201 #endif
202
203 template<typename Derived>
204 bool solveInPlace(MatrixBase<Derived> &bAndX) const;
205
206 LDLT& compute(const MatrixType& matrix);
207
208 template <typename Derived>
209 LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
210
211 /** \returns the internal LDLT decomposition matrix
212 *
213 * TODO: document the storage layout
214 */
215 inline const MatrixType& matrixLDLT() const
216 {
217 eigen_assert(m_isInitialized && "LDLT is not initialized.");
218 return m_matrix;
219 }
220
221 MatrixType reconstructedMatrix() const;
222
223 inline Index rows() const { return m_matrix.rows(); }
224 inline Index cols() const { return m_matrix.cols(); }
225
226 /** \brief Reports whether previous computation was successful.
227 *
228 * \returns \c Success if computation was succesful,
229 * \c NumericalIssue if the matrix.appears to be negative.
230 */
231 ComputationInfo info() const
232 {
233 eigen_assert(m_isInitialized && "LDLT is not initialized.");
234 return Success;
235 }
236
237 protected:
238
239 static void check_template_parameters()
240 {
241 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
242 }
243
244 /** \internal
245 * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
246 * The strict upper part is used during the decomposition, the strict lower
247 * part correspond to the coefficients of L (its diagonal is equal to 1 and
248 * is not stored), and the diagonal entries correspond to D.
249 */
250 MatrixType m_matrix;
251 TranspositionType m_transpositions;
252 TmpMatrixType m_temporary;
253 internal::SignMatrix m_sign;
254 bool m_isInitialized;
255};
256
257namespace internal {
258
259template<int UpLo> struct ldlt_inplace;
260
261template<> struct ldlt_inplace<Lower>
262{
263 template<typename MatrixType, typename TranspositionType, typename Workspace>
264 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
265 {
266 using std::abs;
267 typedef typename MatrixType::Scalar Scalar;
268 typedef typename MatrixType::RealScalar RealScalar;
269 typedef typename MatrixType::Index Index;
270 eigen_assert(mat.rows()==mat.cols());
271 const Index size = mat.rows();
272
273 if (size <= 1)
274 {
275 transpositions.setIdentity();
276 if (numext::real(mat.coeff(0,0)) > 0) sign = PositiveSemiDef;
277 else if (numext::real(mat.coeff(0,0)) < 0) sign = NegativeSemiDef;
278 else sign = ZeroSign;
279 return true;
280 }
281
282 for (Index k = 0; k < size; ++k)
283 {
284 // Find largest diagonal element
285 Index index_of_biggest_in_corner;
286 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
287 index_of_biggest_in_corner += k;
288
289 transpositions.coeffRef(k) = index_of_biggest_in_corner;
290 if(k != index_of_biggest_in_corner)
291 {
292 // apply the transposition while taking care to consider only
293 // the lower triangular part
294 Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
295 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
296 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
297 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
298 for(int i=k+1;i<index_of_biggest_in_corner;++i)
299 {
300 Scalar tmp = mat.coeffRef(i,k);
301 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
302 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
303 }
304 if(NumTraits<Scalar>::IsComplex)
305 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
306 }
307
308 // partition the matrix:
309 // A00 | - | -
310 // lu = A10 | A11 | -
311 // A20 | A21 | A22
312 Index rs = size - k - 1;
313 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
314 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
315 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
316
317 if(k>0)
318 {
319 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
320 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
321 if(rs>0)
322 A21.noalias() -= A20 * temp.head(k);
323 }
324
325 // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
326 // was smaller than the cutoff value. However, soince LDLT is not rank-revealing
327 // we should only make sure we do not introduce INF or NaN values.
328 // LAPACK also uses 0 as the cutoff value.
329 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
330 if((rs>0) && (abs(realAkk) > RealScalar(0)))
331 A21 /= realAkk;
332
333 if (sign == PositiveSemiDef) {
334 if (realAkk < 0) sign = Indefinite;
335 } else if (sign == NegativeSemiDef) {
336 if (realAkk > 0) sign = Indefinite;
337 } else if (sign == ZeroSign) {
338 if (realAkk > 0) sign = PositiveSemiDef;
339 else if (realAkk < 0) sign = NegativeSemiDef;
340 }
341 }
342
343 return true;
344 }
345
346 // Reference for the algorithm: Davis and Hager, "Multiple Rank
347 // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
348 // Trivial rearrangements of their computations (Timothy E. Holy)
349 // allow their algorithm to work for rank-1 updates even if the
350 // original matrix is not of full rank.
351 // Here only rank-1 updates are implemented, to reduce the
352 // requirement for intermediate storage and improve accuracy
353 template<typename MatrixType, typename WDerived>
354 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
355 {
356 using numext::isfinite;
357 typedef typename MatrixType::Scalar Scalar;
358 typedef typename MatrixType::RealScalar RealScalar;
359 typedef typename MatrixType::Index Index;
360
361 const Index size = mat.rows();
362 eigen_assert(mat.cols() == size && w.size()==size);
363
364 RealScalar alpha = 1;
365
366 // Apply the update
367 for (Index j = 0; j < size; j++)
368 {
369 // Check for termination due to an original decomposition of low-rank
370 if (!(isfinite)(alpha))
371 break;
372
373 // Update the diagonal terms
374 RealScalar dj = numext::real(mat.coeff(j,j));
375 Scalar wj = w.coeff(j);
376 RealScalar swj2 = sigma*numext::abs2(wj);
377 RealScalar gamma = dj*alpha + swj2;
378
379 mat.coeffRef(j,j) += swj2/alpha;
380 alpha += swj2/dj;
381
382
383 // Update the terms of L
384 Index rs = size-j-1;
385 w.tail(rs) -= wj * mat.col(j).tail(rs);
386 if(gamma != 0)
387 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
388 }
389 return true;
390 }
391
392 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
393 static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
394 {
395 // Apply the permutation to the input w
396 tmp = transpositions * w;
397
398 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
399 }
400};
401
402template<> struct ldlt_inplace<Upper>
403{
404 template<typename MatrixType, typename TranspositionType, typename Workspace>
405 static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
406 {
407 Transpose<MatrixType> matt(mat);
408 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
409 }
410
411 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
412 static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
413 {
414 Transpose<MatrixType> matt(mat);
415 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
416 }
417};
418
419template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
420{
421 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
422 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
423 static inline MatrixL getL(const MatrixType& m) { return m; }
424 static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
425};
426
427template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
428{
429 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
430 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
431 static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
432 static inline MatrixU getU(const MatrixType& m) { return m; }
433};
434
435} // end namespace internal
436
437/** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
438 */
439template<typename MatrixType, int _UpLo>
440LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const MatrixType& a)
441{
442 check_template_parameters();
443
444 eigen_assert(a.rows()==a.cols());
445 const Index size = a.rows();
446
447 m_matrix = a;
448
449 m_transpositions.resize(size);
450 m_isInitialized = false;
451 m_temporary.resize(size);
452 m_sign = internal::ZeroSign;
453
454 internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign);
455
456 m_isInitialized = true;
457 return *this;
458}
459
460/** Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
461 * \param w a vector to be incorporated into the decomposition.
462 * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
463 * \sa setZero()
464 */
465template<typename MatrixType, int _UpLo>
466template<typename Derived>
467LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
468{
469 const Index size = w.rows();
470 if (m_isInitialized)
471 {
472 eigen_assert(m_matrix.rows()==size);
473 }
474 else
475 {
476 m_matrix.resize(size,size);
477 m_matrix.setZero();
478 m_transpositions.resize(size);
479 for (Index i = 0; i < size; i++)
480 m_transpositions.coeffRef(i) = i;
481 m_temporary.resize(size);
482 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
483 m_isInitialized = true;
484 }
485
486 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
487
488 return *this;
489}
490
491namespace internal {
492template<typename _MatrixType, int _UpLo, typename Rhs>
493struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
494 : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
495{
496 typedef LDLT<_MatrixType,_UpLo> LDLTType;
497 EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
498
499 template<typename Dest> void evalTo(Dest& dst) const
500 {
501 eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
502 // dst = P b
503 dst = dec().transpositionsP() * rhs();
504
505 // dst = L^-1 (P b)
506 dec().matrixL().solveInPlace(dst);
507
508 // dst = D^-1 (L^-1 P b)
509 // more precisely, use pseudo-inverse of D (see bug 241)
510 using std::abs;
511 using std::max;
512 typedef typename LDLTType::MatrixType MatrixType;
513 typedef typename LDLTType::RealScalar RealScalar;
514 const typename Diagonal<const MatrixType>::RealReturnType vectorD(dec().vectorD());
515 // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
516 // as motivated by LAPACK's xGELSS:
517 // RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() *NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
518 // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
519 // diagonal element is not well justified and to numerical issues in some cases.
520 // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
521 RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
522
523 for (Index i = 0; i < vectorD.size(); ++i) {
524 if(abs(vectorD(i)) > tolerance)
525 dst.row(i) /= vectorD(i);
526 else
527 dst.row(i).setZero();
528 }
529
530 // dst = L^-T (D^-1 L^-1 P b)
531 dec().matrixU().solveInPlace(dst);
532
533 // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
534 dst = dec().transpositionsP().transpose() * dst;
535 }
536};
537}
538
539/** \internal use x = ldlt_object.solve(x);
540 *
541 * This is the \em in-place version of solve().
542 *
543 * \param bAndX represents both the right-hand side matrix b and result x.
544 *
545 * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
546 *
547 * This version avoids a copy when the right hand side matrix b is not
548 * needed anymore.
549 *
550 * \sa LDLT::solve(), MatrixBase::ldlt()
551 */
552template<typename MatrixType,int _UpLo>
553template<typename Derived>
554bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
555{
556 eigen_assert(m_isInitialized && "LDLT is not initialized.");
557 eigen_assert(m_matrix.rows() == bAndX.rows());
558
559 bAndX = this->solve(bAndX);
560
561 return true;
562}
563
564/** \returns the matrix represented by the decomposition,
565 * i.e., it returns the product: P^T L D L^* P.
566 * This function is provided for debug purpose. */
567template<typename MatrixType, int _UpLo>
568MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
569{
570 eigen_assert(m_isInitialized && "LDLT is not initialized.");
571 const Index size = m_matrix.rows();
572 MatrixType res(size,size);
573
574 // P
575 res.setIdentity();
576 res = transpositionsP() * res;
577 // L^* P
578 res = matrixU() * res;
579 // D(L^*P)
580 res = vectorD().real().asDiagonal() * res;
581 // L(DL^*P)
582 res = matrixL() * res;
583 // P^T (LDL^*P)
584 res = transpositionsP().transpose() * res;
585
586 return res;
587}
588
589/** \cholesky_module
590 * \returns the Cholesky decomposition with full pivoting without square root of \c *this
591 */
592template<typename MatrixType, unsigned int UpLo>
593inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
594SelfAdjointView<MatrixType, UpLo>::ldlt() const
595{
596 return LDLT<PlainObject,UpLo>(m_matrix);
597}
598
599/** \cholesky_module
600 * \returns the Cholesky decomposition with full pivoting without square root of \c *this
601 */
602template<typename Derived>
603inline const LDLT<typename MatrixBase<Derived>::PlainObject>
604MatrixBase<Derived>::ldlt() const
605{
606 return LDLT<PlainObject>(derived());
607}
608
609} // end namespace Eigen
610
611#endif // EIGEN_LDLT_H
Note: See TracBrowser for help on using the repository browser.