Public Member Functions | Private Types | Private Attributes | List of all members
Spectra::DenseCholesky< Scalar, Uplo > Class Template Reference

#include <DenseCholesky.h>

Public Member Functions

Index cols () const
 
 DenseCholesky (ConstGenericMatrix &mat)
 
int info () const
 
void lower_triangular_solve (const Scalar *x_in, Scalar *y_out) const
 
Index rows () const
 
void upper_triangular_solve (const Scalar *x_in, Scalar *y_out) const
 

Private Types

typedef Eigen::Index Index
 
typedef Eigen::Map< const MatrixMapConstMat
 
typedef Eigen::Map< const VectorMapConstVec
 
typedef Eigen::Map< VectorMapVec
 
typedef Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::DynamicMatrix
 
typedef Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
 

Private Attributes

const typedef Eigen::Ref< const MatrixConstGenericMatrix
 
Eigen::LLT< Matrix, Uplo > m_decomp
 
int m_info
 
const Index m_n
 

Detailed Description

template<typename Scalar, int Uplo = Eigen::Lower>
class Spectra::DenseCholesky< Scalar, Uplo >

This class defines the operations related to Cholesky decomposition on a positive definite matrix, $B=LL'$, where $L$ is a lower triangular matrix. It is mainly used in the SymGEigsSolver generalized eigen solver in the Cholesky decomposition mode.

Definition at line 26 of file DenseCholesky.h.

Member Typedef Documentation

◆ Index

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Index Spectra::DenseCholesky< Scalar, Uplo >::Index
private

Definition at line 29 of file DenseCholesky.h.

◆ MapConstMat

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Map<const Matrix> Spectra::DenseCholesky< Scalar, Uplo >::MapConstMat
private

Definition at line 32 of file DenseCholesky.h.

◆ MapConstVec

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Map<const Vector> Spectra::DenseCholesky< Scalar, Uplo >::MapConstVec
private

Definition at line 33 of file DenseCholesky.h.

◆ MapVec

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Map<Vector> Spectra::DenseCholesky< Scalar, Uplo >::MapVec
private

Definition at line 34 of file DenseCholesky.h.

◆ Matrix

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Spectra::DenseCholesky< Scalar, Uplo >::Matrix
private

Definition at line 30 of file DenseCholesky.h.

◆ Vector

template<typename Scalar , int Uplo = Eigen::Lower>
typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> Spectra::DenseCholesky< Scalar, Uplo >::Vector
private

Definition at line 31 of file DenseCholesky.h.

Constructor & Destructor Documentation

◆ DenseCholesky()

template<typename Scalar , int Uplo = Eigen::Lower>
Spectra::DenseCholesky< Scalar, Uplo >::DenseCholesky ( ConstGenericMatrix mat)
inline

Constructor to create the matrix operation object.

Parameters
matAn Eigen matrix object, whose type can be Eigen::Matrix<Scalar, ...> (e.g. Eigen::MatrixXd and Eigen::MatrixXf), or its mapped version (e.g. Eigen::Map<Eigen::MatrixXd>).

Definition at line 50 of file DenseCholesky.h.

Member Function Documentation

◆ cols()

template<typename Scalar , int Uplo = Eigen::Lower>
Index Spectra::DenseCholesky< Scalar, Uplo >::cols ( ) const
inline

Returns the number of columns of the underlying matrix.

Definition at line 69 of file DenseCholesky.h.

◆ info()

template<typename Scalar , int Uplo = Eigen::Lower>
int Spectra::DenseCholesky< Scalar, Uplo >::info ( ) const
inline

Returns the status of the computation. The full list of enumeration values can be found in Enumerations.

Definition at line 75 of file DenseCholesky.h.

◆ lower_triangular_solve()

template<typename Scalar , int Uplo = Eigen::Lower>
void Spectra::DenseCholesky< Scalar, Uplo >::lower_triangular_solve ( const Scalar x_in,
Scalar y_out 
) const
inline

Performs the lower triangular solving operation $y=L^{-1}x$.

Parameters
x_inPointer to the $x$ vector.
y_outPointer to the $y$ vector.

Definition at line 84 of file DenseCholesky.h.

◆ rows()

template<typename Scalar , int Uplo = Eigen::Lower>
Index Spectra::DenseCholesky< Scalar, Uplo >::rows ( ) const
inline

Returns the number of rows of the underlying matrix.

Definition at line 65 of file DenseCholesky.h.

◆ upper_triangular_solve()

template<typename Scalar , int Uplo = Eigen::Lower>
void Spectra::DenseCholesky< Scalar, Uplo >::upper_triangular_solve ( const Scalar x_in,
Scalar y_out 
) const
inline

Performs the upper triangular solving operation $y=(L')^{-1}x$.

Parameters
x_inPointer to the $x$ vector.
y_outPointer to the $y$ vector.

Definition at line 98 of file DenseCholesky.h.

Member Data Documentation

◆ ConstGenericMatrix

template<typename Scalar , int Uplo = Eigen::Lower>
const typedef Eigen::Ref<const Matrix> Spectra::DenseCholesky< Scalar, Uplo >::ConstGenericMatrix
private

Definition at line 35 of file DenseCholesky.h.

◆ m_decomp

template<typename Scalar , int Uplo = Eigen::Lower>
Eigen::LLT<Matrix, Uplo> Spectra::DenseCholesky< Scalar, Uplo >::m_decomp
private

Definition at line 38 of file DenseCholesky.h.

◆ m_info

template<typename Scalar , int Uplo = Eigen::Lower>
int Spectra::DenseCholesky< Scalar, Uplo >::m_info
private

Definition at line 39 of file DenseCholesky.h.

◆ m_n

template<typename Scalar , int Uplo = Eigen::Lower>
const Index Spectra::DenseCholesky< Scalar, Uplo >::m_n
private

Definition at line 37 of file DenseCholesky.h.


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


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:16:52