Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex > Class Template Reference

#include <Spectra/MatOp/SparseCholesky.h>

Public Types

using Scalar = Scalar_
 

Public Member Functions

template<typename Derived >
 SparseCholesky (const Eigen::SparseMatrixBase< Derived > &mat)
 
Index rows () const
 
Index cols () const
 
CompInfo info () const
 
void lower_triangular_solve (const Scalar *x_in, Scalar *y_out) const
 
void upper_triangular_solve (const Scalar *x_in, Scalar *y_out) const
 

Detailed Description

template<typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
class Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >

This class defines the operations related to Cholesky decomposition on a sparse 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.

Template Parameters
Scalar_The element type of the matrix, for example, float, double, and long double.
UploEither Eigen::Lower or Eigen::Upper, indicating which triangular part of the matrix is used.
FlagsEither Eigen::ColMajor or Eigen::RowMajor, indicating the storage format of the input matrix.
StorageIndexThe type of the indices for the sparse matrix.

Definition at line 36 of file SparseCholesky.h.

Member Typedef Documentation

◆ Scalar

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
using Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::Scalar = Scalar_

Element type of the matrix.

Definition at line 42 of file SparseCholesky.h.

Constructor & Destructor Documentation

◆ SparseCholesky()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
template<typename Derived >
Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::SparseCholesky ( const Eigen::SparseMatrixBase< Derived > &  mat)
inline

Constructor to create the matrix operation object.

Parameters
matAn Eigen sparse matrix object, whose type can be Eigen::SparseMatrix<Scalar, ...> or its mapped version Eigen::Map<Eigen::SparseMatrix<Scalar, ...> >.

Definition at line 64 of file SparseCholesky.h.

Member Function Documentation

◆ rows()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
Index Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::rows ( ) const
inline

Returns the number of rows of the underlying matrix.

Definition at line 83 of file SparseCholesky.h.

◆ cols()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
Index Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::cols ( ) const
inline

Returns the number of columns of the underlying matrix.

Definition at line 87 of file SparseCholesky.h.

◆ info()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
CompInfo Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::info ( ) const
inline

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

Definition at line 93 of file SparseCholesky.h.

◆ lower_triangular_solve()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
void Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::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 102 of file SparseCholesky.h.

◆ upper_triangular_solve()

template<typename Scalar_ , int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
void Spectra::SparseCholesky< Scalar_, Uplo, Flags, StorageIndex >::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 117 of file SparseCholesky.h.


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