Go to the documentation of this file.
19 template<
typename _MatrixType,
int _UpLo>
struct traits<
LDLT<_MatrixType, _UpLo> >
59 template<
typename _MatrixType,
int _UpLo>
class LDLT
112 template<
typename InputType>
129 template<
typename InputType>
149 inline typename Traits::MatrixU
matrixU()
const
156 inline typename Traits::MatrixL
matrixL()
const
191 #ifdef EIGEN_PARSED_BY_DOXYGEN
207 template<
typename Rhs>
212 template<
typename Derived>
215 template<
typename InputType>
227 template <
typename Derived>
263 #ifndef EIGEN_PARSED_BY_DOXYGEN
264 template<
typename RhsType,
typename DstType>
265 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
267 template<
bool Conjugate,
typename RhsType,
typename DstType>
299 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
305 typedef typename TranspositionType::StorageIndex IndexType;
308 bool found_zero_pivot =
false;
313 transpositions.setIdentity();
324 Index index_of_biggest_in_corner;
325 mat.diagonal().tail(
size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
326 index_of_biggest_in_corner += k;
328 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
329 if(k != index_of_biggest_in_corner)
334 mat.row(k).head(k).swap(
mat.row(index_of_biggest_in_corner).head(k));
335 mat.col(k).tail(
s).swap(
mat.col(index_of_biggest_in_corner).tail(
s));
336 std::swap(
mat.coeffRef(k,k),
mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
337 for(
Index i=k+1;
i<index_of_biggest_in_corner;++
i)
344 mat.coeffRef(index_of_biggest_in_corner,k) =
numext::conj(
mat.coeff(index_of_biggest_in_corner,k));
358 temp.head(k) =
mat.diagonal().real().head(k).asDiagonal() *
A10.adjoint();
361 A21.noalias() -= A20 * temp.head(k);
371 if(k==0 && !pivot_is_valid)
378 transpositions.coeffRef(
j) = IndexType(
j);
384 if((rs>0) && pivot_is_valid)
389 if(found_zero_pivot && pivot_is_valid)
ret =
false;
390 else if(!pivot_is_valid) found_zero_pivot =
true;
412 template<
typename MatrixType,
typename WDerived>
443 w.tail(rs) -= wj *
mat.col(
j).tail(rs);
450 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
454 tmp = transpositions *
w;
462 template<
typename MatrixType,
typename TranspositionType,
typename Workspace>
469 template<
typename MatrixType,
typename TranspositionType,
typename Workspace,
typename WType>
497 template<
typename MatrixType,
int _UpLo>
498 template<
typename InputType>
501 check_template_parameters();
506 m_matrix =
a.derived();
514 abs_col_sum = m_matrix.col(
col).tail(
size -
col).template lpNorm<1>() + m_matrix.row(
col).head(
col).template lpNorm<1>();
516 abs_col_sum = m_matrix.col(
col).head(
col).template lpNorm<1>() + m_matrix.row(
col).tail(
size -
col).template lpNorm<1>();
517 if (abs_col_sum > m_l1_norm)
518 m_l1_norm = abs_col_sum;
521 m_transpositions.resize(
size);
522 m_isInitialized =
false;
523 m_temporary.resize(
size);
528 m_isInitialized =
true;
537 template<
typename MatrixType,
int _UpLo>
538 template<
typename Derived>
551 m_transpositions.resize(
size);
553 m_transpositions.coeffRef(
i) = IndexType(
i);
554 m_temporary.resize(
size);
556 m_isInitialized =
true;
564 #ifndef EIGEN_PARSED_BY_DOXYGEN
565 template<
typename _MatrixType,
int _UpLo>
566 template<
typename RhsType,
typename DstType>
569 _solve_impl_transposed<true>(rhs, dst);
572 template<
typename _MatrixType,
int _UpLo>
573 template<
bool Conjugate,
typename RhsType,
typename DstType>
577 dst = m_transpositions * rhs;
581 matrixL().template conjugateIf<!Conjugate>().solveInPlace(dst);
596 for (
Index i = 0;
i < vecD.size(); ++
i)
598 if(
abs(vecD(
i)) > tolerance)
599 dst.row(
i) /= vecD(
i);
601 dst.row(
i).setZero();
606 matrixL().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
610 dst = m_transpositions.transpose() * dst;
627 template<
typename MatrixType,
int _UpLo>
628 template<
typename Derived>
631 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
634 bAndX = this->solve(bAndX);
642 template<
typename MatrixType,
int _UpLo>
645 eigen_assert(m_isInitialized &&
"LDLT is not initialized.");
651 res = transpositionsP() *
res;
655 res = vectorD().real().asDiagonal() *
res;
659 res = transpositionsP().transpose() *
res;
668 template<
typename MatrixType,
unsigned int UpLo>
679 template<
typename Derived>
688 #endif // EIGEN_LDLT_H
LDLT(Index size)
Default Constructor with memory preallocation.
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
#define EIGEN_DEVICE_FUNC
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Namespace containing all symbols from the Eigen library.
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationType
static MatrixU getU(const MatrixType &m)
internal::traits< LDLT< _MatrixType, _UpLo > >::Scalar Scalar
Expression of a fixed-size or dynamic-size block.
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
SolverStorage StorageKind
Eigen::Index Index
The interface type of indices.
MatrixType reconstructedMatrix() const
TranspositionType m_transpositions
Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE bool() isfinite(const Eigen::bfloat16 &h)
const LDLT< PlainObject, UpLo > ldlt() const
static MatrixL getL(const MatrixType &m)
void _solve_impl(const RhsType &rhs, DstType &dst) const
const typedef TriangularView< const MatrixType, UnitUpper > MatrixU
static void check_template_parameters()
const EIGEN_DEVICE_FUNC SignReturnType sign() const
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
bool isNegative(void) const
Traits::MatrixU matrixU() const
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
Expression of the transpose of a matrix.
static const double sigma
LDLT()
Default Constructor.
LDLT & rankUpdate(const MatrixBase< Derived > &w, const RealScalar &alpha=1)
const MatrixType & matrixLDLT() const
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
LDLT & compute(const EigenBase< InputType > &matrix)
internal::SignMatrix m_sign
Traits::MatrixL matrixL() const
Diagonal< const MatrixType > vectorD() const
static const Eigen::internal::all_t all
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
const LDLT & adjoint() const
#define EIGEN_STRONG_INLINE
IndicesType::Scalar StorageIndex
void swap(GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &a, GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &b)
EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & derived()
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
Robust Cholesky decomposition of a matrix with pivoting.
static bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
AnnoyingScalar conj(const AnnoyingScalar &x)
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
const typedef TriangularView< const typename MatrixType::AdjointReturnType, UnitLower > MatrixL
NumTraits< Scalar >::Real RealScalar
static EIGEN_STRONG_INLINE bool update(MatrixType &mat, TranspositionType &transpositions, Workspace &tmp, WType &w, const typename MatrixType::RealScalar &sigma=1)
static const double A10[]
const typedef TriangularView< const typename MatrixType::AdjointReturnType, UnitUpper > MatrixU
EIGEN_DEVICE_FUNC bool abs2(bool x)
Pseudo expression representing a solving operation.
bool solveInPlace(MatrixBase< Derived > &bAndX) const
TmpMatrixType m_temporary
static EIGEN_STRONG_INLINE bool unblocked(MatrixType &mat, TranspositionType &transpositions, Workspace &temp, SignMatrix &sign)
static bool update(MatrixType &mat, const TranspositionType &transpositions, Workspace &tmp, const WType &w, const typename MatrixType::RealScalar &sigma=1)
internal::LDLT_Traits< MatrixType, UpLo > Traits
static MatrixL getL(const MatrixType &m)
static MatrixU getU(const MatrixType &m)
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Base class for all dense matrices, vectors, and expressions.
EIGEN_CONSTEXPR Index size(const T &x)
ComputationInfo info() const
Reports whether previous computation was successful.
const typedef TriangularView< const MatrixType, UnitLower > MatrixL
static bool updateInPlace(MatrixType &mat, MatrixBase< WDerived > &w, const typename MatrixType::RealScalar &sigma=1)
const LDLT< PlainObject > ldlt() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
const Solve< LDLT< _MatrixType, _UpLo >, Rhs > solve(const MatrixBase< Rhs > &b) const
const TranspositionType & transpositionsP() const
A base class for matrix decomposition and solvers.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:11:54