16 template<
typename MatrixType,
int UpLo>
struct LLT_Traits;
50 template<
typename _MatrixType,
int _UpLo>
class LLT 53 typedef _MatrixType MatrixType;
55 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57 Options = MatrixType::Options,
58 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
60 typedef typename MatrixType::Scalar Scalar;
62 typedef typename MatrixType::Index Index;
65 PacketSize = internal::packet_traits<Scalar>::size,
66 AlignmentMask = int(PacketSize)-1,
70 typedef internal::LLT_Traits<MatrixType,UpLo> Traits;
78 LLT() : m_matrix(), m_isInitialized(false) {}
86 LLT(Index size) : m_matrix(size, size),
87 m_isInitialized(false) {}
89 LLT(
const MatrixType& matrix)
90 : m_matrix(matrix.rows(), matrix.cols()),
91 m_isInitialized(
false)
97 inline typename Traits::MatrixU
matrixU()
const 99 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
100 return Traits::getU(m_matrix);
104 inline typename Traits::MatrixL
matrixL()
const 106 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
107 return Traits::getL(m_matrix);
120 template<
typename Rhs>
121 inline const internal::solve_retval<LLT, Rhs>
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());
130 #ifdef EIGEN2_SUPPORT 131 template<
typename OtherDerived,
typename ResultType>
134 *result = this->solve(b);
138 bool isPositiveDefinite()
const {
return true; }
141 template<
typename Derived>
144 LLT& compute(
const MatrixType& matrix);
152 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
156 MatrixType reconstructedMatrix()
const;
166 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
170 inline Index rows()
const {
return m_matrix.rows(); }
171 inline Index cols()
const {
return m_matrix.cols(); }
173 template<
typename VectorType>
174 LLT rankUpdate(
const VectorType& vec,
const RealScalar& sigma = 1);
182 bool m_isInitialized;
188 template<
typename Scalar,
int UpLo>
struct llt_inplace;
190 template<
typename MatrixType,
typename VectorType>
191 static typename MatrixType::Index llt_rank_update_lower(MatrixType& mat,
const VectorType& vec,
const typename MatrixType::RealScalar& sigma)
194 typedef typename MatrixType::Scalar Scalar;
195 typedef typename MatrixType::RealScalar RealScalar;
196 typedef typename MatrixType::Index Index;
197 typedef typename MatrixType::ColXpr ColXpr;
198 typedef typename internal::remove_all<ColXpr>::type ColXprCleaned;
199 typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
201 typedef typename TempVectorType::SegmentReturnType TempVecSegment;
203 Index n = mat.cols();
204 eigen_assert(mat.rows()==n && vec.size()==n);
213 temp = sqrt(sigma) * vec;
215 for(Index i=0; i<n; ++i)
223 ColXprSegment x(mat.col(i).tail(rs));
224 TempVecSegment y(temp.tail(rs));
225 apply_rotation_in_the_plane(x, y, g);
233 for(Index j=0; j<n; ++j)
235 RealScalar Ljj = numext::real(mat.coeff(j,j));
236 RealScalar dj = numext::abs2(Ljj);
237 Scalar wj = temp.coeff(j);
238 RealScalar swj2 = sigma*numext::abs2(wj);
239 RealScalar gamma = dj*beta + swj2;
241 RealScalar x = dj + swj2/beta;
242 if (x<=RealScalar(0))
244 RealScalar nLjj = sqrt(x);
245 mat.coeffRef(j,j) = nLjj;
252 temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
254 mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
261 template<
typename Scalar>
struct llt_inplace<Scalar, Lower>
264 template<
typename MatrixType>
265 static typename MatrixType::Index unblocked(MatrixType& mat)
268 typedef typename MatrixType::Index Index;
270 eigen_assert(mat.rows()==mat.cols());
271 const Index size = mat.rows();
272 for(Index k = 0; k < size; ++k)
280 RealScalar x = numext::real(mat.coeff(k,k));
281 if (k>0) x -= A10.squaredNorm();
282 if (x<=RealScalar(0))
284 mat.coeffRef(k,k) = x = sqrt(x);
285 if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
286 if (rs>0) A21 *= RealScalar(1)/x;
291 template<
typename MatrixType>
292 static typename MatrixType::Index blocked(MatrixType& m)
294 typedef typename MatrixType::Index Index;
295 eigen_assert(m.rows()==m.cols());
296 Index size = m.rows();
300 Index blockSize = size/8;
301 blockSize = (blockSize/16)*16;
302 blockSize = (std::min)((std::max)(blockSize,Index(8)), Index(128));
304 for (Index k=0; k<size; k+=blockSize)
310 Index bs = (std::min)(blockSize, size-k);
311 Index rs = size - k - bs;
317 if((ret=unblocked(A11))>=0)
return k+ret;
318 if(rs>0) A11.adjoint().template triangularView<Upper>().
template solveInPlace<OnTheRight>(A21);
319 if(rs>0) A22.template selfadjointView<Lower>().rankUpdate(A21,-1);
324 template<
typename MatrixType,
typename VectorType>
325 static typename MatrixType::Index rankUpdate(MatrixType& mat,
const VectorType& vec,
const RealScalar& sigma)
327 return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
331 template<
typename Scalar>
struct llt_inplace<Scalar, Upper>
335 template<
typename MatrixType>
336 static EIGEN_STRONG_INLINE
typename MatrixType::Index unblocked(MatrixType& mat)
339 return llt_inplace<Scalar, Lower>::unblocked(matt);
341 template<
typename MatrixType>
342 static EIGEN_STRONG_INLINE
typename MatrixType::Index blocked(MatrixType& mat)
345 return llt_inplace<Scalar, Lower>::blocked(matt);
347 template<
typename MatrixType,
typename VectorType>
348 static typename MatrixType::Index rankUpdate(MatrixType& mat,
const VectorType& vec,
const RealScalar& sigma)
351 return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
355 template<
typename MatrixType>
struct LLT_Traits<MatrixType,Lower>
359 static inline MatrixL getL(
const MatrixType& m) {
return m; }
360 static inline MatrixU getU(
const MatrixType& m) {
return m.
adjoint(); }
361 static bool inplace_decomposition(MatrixType& m)
362 {
return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m)==-1; }
365 template<
typename MatrixType>
struct LLT_Traits<MatrixType,Upper>
369 static inline MatrixL getL(
const MatrixType& m) {
return m.
adjoint(); }
370 static inline MatrixU getU(
const MatrixType& m) {
return m; }
371 static bool inplace_decomposition(MatrixType& m)
372 {
return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m)==-1; }
384 template<
typename MatrixType,
int _UpLo>
387 eigen_assert(a.rows()==a.cols());
388 const Index size = a.rows();
389 m_matrix.resize(size, size);
392 m_isInitialized =
true;
393 bool ok = Traits::inplace_decomposition(m_matrix);
404 template<
typename _MatrixType,
int _UpLo>
405 template<
typename VectorType>
408 EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
409 eigen_assert(v.size()==m_matrix.cols());
410 eigen_assert(m_isInitialized);
411 if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
420 template<
typename _MatrixType,
int UpLo,
typename Rhs>
421 struct solve_retval<LLT<_MatrixType, UpLo>, Rhs>
422 : solve_retval_base<LLT<_MatrixType, UpLo>, Rhs>
425 EIGEN_MAKE_SOLVE_HELPERS(LLTType,Rhs)
427 template<
typename Dest>
void evalTo(Dest& dst)
const 430 dec().solveInPlace(dst);
448 template<
typename MatrixType,
int _UpLo>
449 template<
typename Derived>
452 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
453 eigen_assert(m_matrix.rows()==bAndX.rows());
454 matrixL().solveInPlace(bAndX);
455 matrixU().solveInPlace(bAndX);
461 template<
typename MatrixType,
int _UpLo>
464 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
465 return matrixL() * matrixL().adjoint().toDenseMatrix();
471 template<
typename Derived>
481 template<
typename MatrixType,
unsigned int UpLo>
490 #endif // EIGEN_LLT_H const internal::solve_retval< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LLT.h:122
LLT()
Default Constructor.
Definition: LLT.h:78
Expression of the transpose of a matrix.
Definition: Transpose.h:57
Definition: Constants.h:378
const TriangularView< const typename MatrixType::AdjointReturnType, TransposeMode > adjoint() const
Definition: TriangularMatrix.h:264
LLT(Index size)
Default Constructor with memory preallocation.
Definition: LLT.h:86
Rotation given by a cosine-sine pair.
Definition: ForwardDeclarations.h:228
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const LLT< PlainObject > llt() const
Definition: LLT.h:473
MatrixType reconstructedMatrix() const
Definition: LLT.h:462
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:483
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Definition: Jacobi.h:148
LLT & compute(const MatrixType &matrix)
Definition: LLT.h:385
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:50
const MatrixType & matrixLLT() const
Definition: LLT.h:150
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LLT.h:164
Definition: Eigen_Colamd.h:54
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:102
Base class for triangular part in a matrix.
Definition: TriangularMatrix.h:158
Definition: Constants.h:376
ComputationInfo
Definition: Constants.h:374
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Traits::MatrixL matrixL() const
Definition: LLT.h:104
Traits::MatrixU matrixU() const
Definition: LLT.h:97