16 template<
typename MatrixType,
int UpLo>
struct LLT_Traits;
56 template<
typename _MatrixType,
int _UpLo>
class LLT 61 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
62 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
63 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
72 AlignmentMask =
int(PacketSize)-1,
84 LLT() : m_matrix(), m_isInitialized(false) {}
92 explicit LLT(Index
size) : m_matrix(size, size),
93 m_isInitialized(false) {}
95 template<
typename InputType>
97 : m_matrix(matrix.rows(), matrix.cols()),
98 m_isInitialized(false)
110 template<
typename InputType>
112 : m_matrix(matrix.derived()),
113 m_isInitialized(false)
119 inline typename Traits::MatrixU
matrixU()
const 121 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
122 return Traits::getU(m_matrix);
126 inline typename Traits::MatrixL
matrixL()
const 128 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
129 return Traits::getL(m_matrix);
142 template<
typename Rhs>
146 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
148 &&
"LLT::solve(): invalid number of rows of the right hand side matrix b");
152 template<
typename Derived>
155 template<
typename InputType>
163 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
174 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
178 MatrixType reconstructedMatrix()
const;
188 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
199 inline Index
rows()
const {
return m_matrix.rows(); }
200 inline Index
cols()
const {
return m_matrix.cols(); }
202 template<
typename VectorType>
203 LLT rankUpdate(
const VectorType& vec,
const RealScalar& sigma = 1);
205 #ifndef EIGEN_PARSED_BY_DOXYGEN 206 template<
typename RhsType,
typename DstType>
208 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
232 template<
typename MatrixType,
typename VectorType>
254 temp =
sqrt(sigma) * vec;
264 ColXprSegment
x(mat.col(i).tail(rs));
265 TempVecSegment
y(temp.tail(rs));
278 Scalar wj = temp.coeff(j);
280 RealScalar gamma = dj*beta + swj2;
282 RealScalar
x = dj + swj2/beta;
285 RealScalar nLjj =
sqrt(x);
286 mat.coeffRef(j,j) = nLjj;
293 temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
295 mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*
numext::conj(wj)/gamma)*temp.tail(rs);
305 template<
typename MatrixType>
321 if (k>0) x -= A10.squaredNorm();
324 mat.coeffRef(k,k) = x =
sqrt(x);
325 if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
331 template<
typename MatrixType>
339 Index blockSize = size/8;
340 blockSize = (blockSize/16)*16;
350 Index rs = size - k - bs;
356 if((ret=unblocked(A11))>=0)
return k+
ret;
357 if(rs>0) A11.adjoint().template triangularView<Upper>().
template solveInPlace<OnTheRight>(A21);
363 template<
typename MatrixType,
typename VectorType>
374 template<
typename MatrixType>
380 template<
typename MatrixType>
386 template<
typename MatrixType,
typename VectorType>
423 template<
typename MatrixType,
int _UpLo>
424 template<
typename InputType>
427 check_template_parameters();
431 m_matrix.resize(size, size);
441 abs_col_sum = m_matrix.col(
col).tail(size -
col).template lpNorm<1>() + m_matrix.row(
col).head(
col).template lpNorm<1>();
443 abs_col_sum = m_matrix.col(
col).head(
col).template lpNorm<1>() + m_matrix.row(
col).tail(size -
col).template lpNorm<1>();
444 if (abs_col_sum > m_l1_norm)
445 m_l1_norm = abs_col_sum;
448 m_isInitialized =
true;
449 bool ok = Traits::inplace_decomposition(m_matrix);
460 template<
typename _MatrixType,
int _UpLo>
461 template<
typename VectorType>
475 #ifndef EIGEN_PARSED_BY_DOXYGEN 476 template<
typename _MatrixType,
int _UpLo>
477 template<
typename RhsType,
typename DstType>
498 template<
typename MatrixType,
int _UpLo>
499 template<
typename Derived>
502 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
504 matrixL().solveInPlace(bAndX);
505 matrixU().solveInPlace(bAndX);
511 template<
typename MatrixType,
int _UpLo>
514 eigen_assert(m_isInitialized &&
"LLT is not initialized.");
515 return matrixL() * matrixL().adjoint().toDenseMatrix();
522 template<
typename Derived>
533 template<
typename MatrixType,
unsigned int UpLo>
542 #endif // EIGEN_LLT_H Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
#define EIGEN_STRONG_INLINE
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
VectorBlock< Derived > SegmentReturnType
MatrixType::StorageIndex StorageIndex
EIGEN_DEVICE_FUNC Index rows() const
const TriangularView< const MatrixType, Lower > MatrixL
bool is_same_dense(const T1 &mat1, const T2 &mat2, typename enable_if< has_direct_access< T1 >::ret &&has_direct_access< T2 >::ret, T1 >::type *=0)
MatrixType reconstructedMatrix() const
Expression of the transpose of a matrix.
EIGEN_DEVICE_FUNC ColXpr col(Index i)
This is the const version of col().
MatrixType::Scalar Scalar
LLT(Index size)
Default Constructor with memory preallocation.
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1, !IsRowMajor > ColXpr
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
static constexpr size_t size(Tuple< Args... > &)
Provides access to the number of elements in a tuple as a compile-time constant expression.
Rotation given by a cosine-sine pair.
static bool inplace_decomposition(MatrixType &m)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
const TriangularView< const typename MatrixType::AdjointReturnType, Upper > MatrixU
LLT(const EigenBase< InputType > &matrix)
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
LLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
static EIGEN_STRONG_INLINE Index blocked(MatrixType &mat)
LLT rankUpdate(const VectorType &vec, const RealScalar &sigma=1)
const LLT & adjoint() const
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
Map< Matrix< Scalar, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > MatrixType
internal::LLT_Traits< MatrixType, UpLo > Traits
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
NumTraits< Scalar >::Real RealScalar
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
static void check_template_parameters()
ComputationInfo info() const
Reports whether previous computation was successful.
static Index rankUpdate(MatrixType &mat, const VectorType &vec, const RealScalar &sigma)
static Index unblocked(MatrixType &mat)
static bool inplace_decomposition(MatrixType &m)
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Traits::MatrixL matrixL() const
static Index llt_rank_update_lower(MatrixType &mat, const VectorType &vec, const typename MatrixType::RealScalar &sigma)
Expression of a fixed-size or dynamic-size block.
static EIGEN_STRONG_INLINE Index unblocked(MatrixType &mat)
static MatrixU getU(const MatrixType &m)
Traits::MatrixU matrixU() const
static MatrixU getU(const MatrixType &m)
static Index rankUpdate(MatrixType &mat, const VectorType &vec, const RealScalar &sigma)
Expression of a triangular part in a matrix.
EIGEN_DEVICE_FUNC Index cols() const
void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
const TriangularView< const MatrixType, Upper > MatrixU
bool ok()
global method to check whether to proceed or cancel the current action
const LLT< PlainObject > llt() const
LLT & compute(const EigenBase< InputType > &matrix)
static Index blocked(MatrixType &m)
const TriangularView< const typename MatrixType::AdjointReturnType, Lower > MatrixL
const MatrixType & matrixLLT() const
Pseudo expression representing a solving operation.
LLT()
Default Constructor.
static MatrixL getL(const MatrixType &m)
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
EIGEN_DEVICE_FUNC Derived & derived()
Base class for all dense matrices, vectors, and expressions.
const LLT< PlainObject, UpLo > llt() const
static MatrixL getL(const MatrixType &m)
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
NumTraits< typename MatrixType::Scalar >::Real RealScalar
void solveInPlace(const MatrixBase< Derived > &bAndX) const
NumTraits< Scalar >::Real RealScalar