Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > Class Template Reference

A direct sparse LDLT Cholesky factorizations without square root, for symmetric non-hermitian matrices. More...

#include <SimplicialCholesky.h>

+ Inheritance diagram for Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >:

Public Types

enum  { UpLo = UpLo_ }
 
typedef MatrixType_ MatrixType
 
typedef SimplicialCholeskyBase< SimplicialNonHermitianLDLTBase
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
typedef internal::traits< SimplicialNonHermitianLDLTTraits
 
typedef Traits::MatrixL MatrixL
 
typedef Traits::MatrixU MatrixU
 
- Public Types inherited from Eigen::SimplicialCholeskyBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
enum  
 
enum  
 
typedef internal::traits< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >::MatrixType MatrixType
 
typedef internal::traits< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >::OrderingType OrderingType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef internal::traits< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >::DiagonalScalar DiagonalScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
typedef CholMatrixType const * ConstCholMatrixPtr
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorI
 

Public Member Functions

 SimplicialNonHermitianLDLT ()
 
 SimplicialNonHermitianLDLT (const MatrixType &matrix)
 
const VectorType vectorD () const
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
SimplicialNonHermitianLDLTcompute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &a)
 
void factorize (const MatrixType &a)
 
Scalar determinant () const
 
- Public Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
 SimplicialCholeskyBase ()
 
 SimplicialCholeskyBase (const MatrixType &matrix)
 
 ~SimplicialCholeskyBase ()
 
SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
Index cols () const
 
Index rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & setShift (const DiagonalScalar &offset, const DiagonalScalar &scale=1)
 
void dumpMemory (Stream &s)
 
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
const Solve< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
void compute (const MatrixType &matrix)
 
void factorize (const MatrixType &a)
 
void factorize_preordered (const CholMatrixType &a)
 
void analyzePattern (const MatrixType &a)
 
void analyzePattern_preordered (const CholMatrixType &a, bool doLDLT)
 
void ordering (const MatrixType &a, ConstCholMatrixPtr &pmat, CholMatrixType &ap)
 
DiagonalScalar getDiag (Scalar x)
 
Scalar getSymm (Scalar x)
 
- Protected Attributes inherited from Eigen::SimplicialCholeskyBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
ComputationInfo m_info
 
bool m_factorizationIsOk
 
bool m_analysisIsOk
 
CholMatrixType m_matrix
 
VectorType m_diag
 
VectorI m_parent
 
VectorI m_workSpace
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_P
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_Pinv
 
DiagonalScalar m_shiftOffset
 
DiagonalScalar m_shiftScale
 
- Protected Attributes inherited from Eigen::SparseSolverBase< SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ > >
bool m_isInitialized
 

Detailed Description

template<typename MatrixType_, int UpLo_, typename Ordering_>
class Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >

A direct sparse LDLT Cholesky factorizations without square root, for symmetric non-hermitian matrices.

This class provides a LDL^T Cholesky factorizations without square root of sparse matrices that are symmetric but not hermitian. For real matrices, this is equivalent to the regular LDLT factorization. The factorization allows for solving A.X = B where X and B can be either dense or sparse.

In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization such that the factorized matrix is P A P^-1.

Template Parameters
MatrixType_the type of the sparse matrix A, it must be a SparseMatrix<>
UpLo_the triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
Ordering_The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>

\implsparsesolverconcept

See also
class SimplicialNonHermitianLLT, SimplicialLDLT, class AMDOrdering, class NaturalOrdering

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef SimplicialCholeskyBase<SimplicialNonHermitianLDLT> Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::Base

◆ CholMatrixType

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef SparseMatrix<Scalar, ColMajor, StorageIndex> Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::CholMatrixType

◆ MatrixL

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Traits::MatrixL Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::MatrixL

◆ MatrixType

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType_ Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::MatrixType

◆ MatrixU

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Traits::MatrixU Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::MatrixU

◆ RealScalar

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::RealScalar Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::RealScalar

◆ Scalar

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::Scalar Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::Scalar

◆ StorageIndex

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::StorageIndex Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::StorageIndex

◆ Traits

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef internal::traits<SimplicialNonHermitianLDLT> Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::Traits

◆ VectorType

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Matrix<Scalar, Dynamic, 1> Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::VectorType

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
anonymous enum
Enumerator
UpLo 
626 { UpLo = UpLo_ };
@ UpLo
Definition: SimplicialCholesky.h:626

Constructor & Destructor Documentation

◆ SimplicialNonHermitianLDLT() [1/2]

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::SimplicialNonHermitianLDLT ( )
inline

Default constructor

639 : Base() {}
SimplicialCholeskyBase< SimplicialNonHermitianLDLT > Base
Definition: SimplicialCholesky.h:627

◆ SimplicialNonHermitianLDLT() [2/2]

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::SimplicialNonHermitianLDLT ( const MatrixType matrix)
inlineexplicit

Constructs and performs the LLT factorization of matrix

642 : Base(matrix) {}
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

Member Function Documentation

◆ analyzePattern()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
void Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::analyzePattern ( const MatrixType a)
inline

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
673 { Base::template analyzePattern<true, true>(a); }
const Scalar * a
Definition: level2_cplx_impl.h:32

References a.

◆ compute()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
SimplicialNonHermitianLDLT& Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::compute ( const MatrixType matrix)
inline

Computes the sparse Cholesky decomposition of matrix

662  {
663  Base::template compute<true, true>(matrix);
664  return *this;
665  }

References matrix().

◆ determinant()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Scalar Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization
684 { return Base::m_diag.prod(); }
VectorType m_diag
Definition: SimplicialCholesky.h:241

References Eigen::SimplicialCholeskyBase< Derived >::m_diag.

◆ factorize()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
void Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::factorize ( const MatrixType a)
inline

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
681 { Base::template factorize<true, true>(a); }

References a.

◆ matrixL()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
const MatrixL Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::matrixL ( ) const
inline
Returns
an expression of the factor L
650  {
651  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
652  return Traits::getL(Base::m_matrix);
653  }
#define eigen_assert(x)
Definition: Macros.h:910
CholMatrixType m_matrix
Definition: SimplicialCholesky.h:240
bool m_factorizationIsOk
Definition: SimplicialCholesky.h:237

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk, and Eigen::SimplicialCholeskyBase< Derived >::m_matrix.

◆ matrixU()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
const MatrixU Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::matrixU ( ) const
inline
Returns
an expression of the factor U (= L^*)
656  {
657  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
658  return Traits::getU(Base::m_matrix);
659  }

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk, and Eigen::SimplicialCholeskyBase< Derived >::m_matrix.

◆ vectorD()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
const VectorType Eigen::SimplicialNonHermitianLDLT< MatrixType_, UpLo_, Ordering_ >::vectorD ( ) const
inline
Returns
a vector expression of the diagonal D
645  {
646  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
647  return Base::m_diag;
648  }

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_diag, and Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk.


The documentation for this class was generated from the following file: