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

#include <UpperHessenbergQR.h>

Inheritance diagram for Spectra::TridiagQR< Scalar >:
Inheritance graph
[legend]

Public Member Functions

void compute (ConstGenericMatrix &mat, const Scalar &shift=Scalar(0)) override
 
void matrix_QtHQ (ComplexMatrix &dest) const
 
void matrix_QtHQ (Matrix &dest) const override
 
Matrix matrix_R () const override
 
 TridiagQR (ConstGenericMatrix &mat, const Scalar &shift=Scalar(0))
 
 TridiagQR (Index size)
 
- Public Member Functions inherited from Spectra::UpperHessenbergQR< double >
void apply_QtY (GenericMatrix Y) const
 
void apply_QtY (Vector &Y) const
 
void apply_QY (GenericMatrix Y) const
 
void apply_QY (Vector &Y) const
 
void apply_YQ (GenericMatrix Y) const
 
void apply_YQt (GenericMatrix Y) const
 
virtual void compute (ConstGenericMatrix &mat, const double &shift=double(0))
 
virtual void matrix_QtHQ (Matrix &dest) const
 
virtual Matrix matrix_R () const
 
 UpperHessenbergQR (ConstGenericMatrix &mat, const double &shift=double(0))
 
 UpperHessenbergQR (Index size)
 
virtual ~UpperHessenbergQR ()
 

Private Types

using ComplexMatrix = Eigen::Matrix< std::complex< Scalar >, Eigen::Dynamic, Eigen::Dynamic >
 
using ConstGenericMatrix = const Eigen::Ref< const Matrix >
 
using Index = Eigen::Index
 
using Matrix = Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic >
 
using Vector = Eigen::Matrix< Scalar, Eigen::Dynamic, 1 >
 

Private Attributes

Vector m_R_diag
 
Vector m_R_supd
 
Vector m_R_supd2
 
Vector m_T_diag
 
Vector m_T_subd
 

Additional Inherited Members

- Static Protected Member Functions inherited from Spectra::UpperHessenbergQR< double >
static void compute_rotation (const double &x, const double &y, double &r, double &c, double &s)
 
static void stable_scaling (const double &a, const double &b, double &r, double &c, double &s)
 
- Protected Attributes inherited from Spectra::UpperHessenbergQR< double >
bool m_computed
 
Index m_n
 
Array m_rot_cos
 
Array m_rot_sin
 
double m_shift
 

Detailed Description

template<typename Scalar = double>
class Spectra::TridiagQR< Scalar >

Perform the QR decomposition of a tridiagonal matrix, a special case of upper Hessenberg matrices.

Template Parameters
ScalarThe element type of the matrix. Currently supported types are float, double and long double.

Definition at line 545 of file UpperHessenbergQR.h.

Member Typedef Documentation

◆ ComplexMatrix

template<typename Scalar = double>
using Spectra::TridiagQR< Scalar >::ComplexMatrix = Eigen::Matrix<std::complex<Scalar>, Eigen::Dynamic, Eigen::Dynamic>
private

Definition at line 552 of file UpperHessenbergQR.h.

◆ ConstGenericMatrix

template<typename Scalar = double>
using Spectra::TridiagQR< Scalar >::ConstGenericMatrix = const Eigen::Ref<const Matrix>
private

Definition at line 551 of file UpperHessenbergQR.h.

◆ Index

template<typename Scalar = double>
using Spectra::TridiagQR< Scalar >::Index = Eigen::Index
private

Definition at line 548 of file UpperHessenbergQR.h.

◆ Matrix

template<typename Scalar = double>
using Spectra::TridiagQR< Scalar >::Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>
private

Definition at line 549 of file UpperHessenbergQR.h.

◆ Vector

template<typename Scalar = double>
using Spectra::TridiagQR< Scalar >::Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>
private

Definition at line 550 of file UpperHessenbergQR.h.

Constructor & Destructor Documentation

◆ TridiagQR() [1/2]

template<typename Scalar = double>
Spectra::TridiagQR< Scalar >::TridiagQR ( Index  size)
inline

Constructor to preallocate memory. Computation can be performed later by calling the compute() method.

Definition at line 571 of file UpperHessenbergQR.h.

◆ TridiagQR() [2/2]

template<typename Scalar = double>
Spectra::TridiagQR< Scalar >::TridiagQR ( ConstGenericMatrix mat,
const Scalar shift = Scalar(0) 
)
inline

Constructor to create an object that performs and stores the QR decomposition of a tridiagonal matrix mat, with an optional shift: $T-sI=QR$. Here $T$ stands for the matrix mat, and $s$ is the shift.

Parameters
matMatrix type can be Eigen::Matrix<Scalar, ...> (e.g. Eigen::MatrixXd and Eigen::MatrixXf), or its mapped version (e.g. Eigen::Map<Eigen::MatrixXd>). Only the diagonal and subdiagonal elements of the matrix are used.

Definition at line 586 of file UpperHessenbergQR.h.

Member Function Documentation

◆ compute()

template<typename Scalar = double>
void Spectra::TridiagQR< Scalar >::compute ( ConstGenericMatrix mat,
const Scalar shift = Scalar(0) 
)
inlineoverride

Compute the QR decomposition of a tridiagonal matrix with an optional shift.

Parameters
matMatrix type can be Eigen::Matrix<Scalar, ...> (e.g. Eigen::MatrixXd and Eigen::MatrixXf), or its mapped version (e.g. Eigen::Map<Eigen::MatrixXd>). Only the diagonal and subdiagonal elements of the matrix are used.

Definition at line 601 of file UpperHessenbergQR.h.

◆ matrix_QtHQ() [1/2]

template<typename Scalar = double>
void Spectra::TridiagQR< Scalar >::matrix_QtHQ ( ComplexMatrix dest) const
inline

The version of matrix_QtHQ() when dest has a complex value type.

This is used in Hermitian eigen solvers where the result is stored as a complex matrix.

Definition at line 787 of file UpperHessenbergQR.h.

◆ matrix_QtHQ() [2/2]

template<typename Scalar = double>
void Spectra::TridiagQR< Scalar >::matrix_QtHQ ( Matrix dest) const
inlineoverride

Overwrite dest with $Q'TQ = RQ + sI$, where $T$ is the input matrix mat, and $s$ is the shift. The result is a tridiagonal matrix.

Parameters
matThe matrix to be overwritten, whose type should be Eigen::Matrix<Scalar, ...>, depending on the template parameter Scalar defined.

Definition at line 713 of file UpperHessenbergQR.h.

◆ matrix_R()

template<typename Scalar = double>
Matrix Spectra::TridiagQR< Scalar >::matrix_R ( ) const
inlineoverride

Return the $R$ matrix in the QR decomposition, which is an upper triangular matrix.

Returns
Returned matrix type will be Eigen::Matrix<Scalar, ...>, depending on the template parameter Scalar defined.

Definition at line 693 of file UpperHessenbergQR.h.

Member Data Documentation

◆ m_R_diag

template<typename Scalar = double>
Vector Spectra::TridiagQR< Scalar >::m_R_diag
private

Definition at line 562 of file UpperHessenbergQR.h.

◆ m_R_supd

template<typename Scalar = double>
Vector Spectra::TridiagQR< Scalar >::m_R_supd
private

Definition at line 563 of file UpperHessenbergQR.h.

◆ m_R_supd2

template<typename Scalar = double>
Vector Spectra::TridiagQR< Scalar >::m_R_supd2
private

Definition at line 564 of file UpperHessenbergQR.h.

◆ m_T_diag

template<typename Scalar = double>
Vector Spectra::TridiagQR< Scalar >::m_T_diag
private

Definition at line 560 of file UpperHessenbergQR.h.

◆ m_T_subd

template<typename Scalar = double>
Vector Spectra::TridiagQR< Scalar >::m_T_subd
private

Definition at line 561 of file UpperHessenbergQR.h.


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


gtsam
Author(s):
autogenerated on Fri Mar 28 2025 03:16:39