TooN
Public Member Functions | List of all members
TooN::SymEigen< Size, Precision > Class Template Reference

Performs eigen decomposition of a matrix. More...

#include <SymEigen.h>

Public Member Functions

 SymEigen (int m)
 Initialise this eigen decomposition but do no immediately perform a decomposition. More...
 
template<int R, int C, typename B >
 SymEigen (const Matrix< R, C, Precision, B > &m)
 Construct the eigen decomposition of a matrix. More...
 
template<int R, int C, typename B >
void compute (const Matrix< R, C, Precision, B > &m)
 Perform the eigen decomposition of a matrix.
 
template<int S, typename P , typename B >
Vector< Size, Precision > backsub (const Vector< S, P, B > &rhs) const
 Calculate result of multiplying the (pseudo-)inverse of M by a vector. More...
 
template<int R, int C, typename P , typename B >
Matrix< Size, C, Precision > backsub (const Matrix< R, C, P, B > &rhs) const
 Calculate result of multiplying the (pseudo-)inverse of M by another matrix. More...
 
Matrix< Size, Size, Precision > get_pinv (const double condition=Internal::symeigen_condition_no) const
 Calculate (pseudo-)inverse of the matrix. More...
 
Vector< Size, Precision > get_inv_diag (const double condition) const
 Calculates the reciprocals of the eigenvalues of the matrix. More...
 
Matrix< Size, Size, Precision > & get_evectors ()
 Returns the eigenvectors of the matrix. More...
 
const Matrix< Size, Size, Precision > & get_evectors () const
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Vector< Size, Precision > & get_evalues ()
 Returns the eigenvalues of the matrix. More...
 
const Vector< Size, Precision > & get_evalues () const
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
bool is_posdef () const
 Is the matrix positive definite?
 
bool is_negdef () const
 Is the matrix negative definite?
 
Precision get_determinant () const
 Get the determinant of the matrix.
 
Matrix< Size, Size, Precision > get_sqrtm () const
 Calculate the square root of a matrix which is a matrix M such that M.T*M=A. More...
 
Matrix< Size, Size, Precision > get_isqrtm (const double condition=Internal::symeigen_condition_no) const
 Calculate the inverse square root of a matrix which is a matrix M such that M.T*M=A^-1. More...
 

Detailed Description

template<int Size = Dynamic, typename Precision = double>
class TooN::SymEigen< Size, Precision >

Performs eigen decomposition of a matrix.

Real symmetric (and hence square matrices) can be decomposed into

\[M = U \times \Lambda \times U^T\]

where \(U\) is an orthogonal matrix (and hence \(U^T = U^{-1}\)) whose columns are the eigenvectors of \(M\) and \(\Lambda\) is a diagonal matrix whose entries are the eigenvalues of \(M\). These quantities are often of use directly, and can be obtained as follows:

// construct M
Matrix<3> M(3,3);
M[0]=makeVector(4,0,2);
M[1]=makeVector(0,5,3);
M[2]=makeVector(2,3,6);
// create the eigen decomposition of M
SymEigen<3> eigM(M);
cout << "A=" << M << endl;
cout << "(E,v)=eig(A)" << endl;
// print the smallest eigenvalue
cout << "v[0]=" << eigM.get_evalues()[0] << endl;
// print the associated eigenvector
cout << "E[0]=" << eigM.get_evectors()[0] << endl;

Further, provided the eigenvalues are nonnegative, the square root of a matrix and its inverse can also be obtained,

// print the square root of the matrix.
cout << "R=sqrtm(A)=" << eigM.get_sqrtm() << endl;
// print the square root of the matrix squared.
cout << "(should equal A), R^T*R="
<< eigM.get_sqrtm().T() * eigM.get_sqrtm() << endl;
// print the inverse of the matrix.
cout << "A^-1=" << eigM.get_pinv() << endl;
// print the inverse square root of the matrix.
cout << "C=isqrtm(A)=" << eigM.get_isqrtm() << endl;
// print the inverse square root of the matrix squared.
cout << "(should equal A^-1), C^T*C="
<< eigM.get_isqrtm().T() * eigM.get_isqrtm() << endl;

This decomposition is very similar to the SVD (q.v.), and can be used to solve equations using backsub() or get_pinv(), with the same treatment of condition numbers.

SymEigen<> (= SymEigen<-1>) can be used to create an eigen decomposition whose size is determined at run-time.

Constructor & Destructor Documentation

◆ SymEigen() [1/2]

template<int Size = Dynamic, typename Precision = double>
TooN::SymEigen< Size, Precision >::SymEigen ( int  m)
inline

Initialise this eigen decomposition but do no immediately perform a decomposition.

Parameters
mThe size of the matrix to perform the eigen decomposition on.

◆ SymEigen() [2/2]

template<int Size = Dynamic, typename Precision = double>
template<int R, int C, typename B >
TooN::SymEigen< Size, Precision >::SymEigen ( const Matrix< R, C, Precision, B > &  m)
inline

Construct the eigen decomposition of a matrix.

This initialises the class, and performs the decomposition immediately.

Member Function Documentation

◆ backsub() [1/2]

template<int Size = Dynamic, typename Precision = double>
template<int S, typename P , typename B >
Vector<Size, Precision> TooN::SymEigen< Size, Precision >::backsub ( const Vector< S, P, B > &  rhs) const
inline

Calculate result of multiplying the (pseudo-)inverse of M by a vector.

For a vector \(b\), this calculates \(M^{\dagger}b\) by back substitution (i.e. without explictly calculating the (pseudo-)inverse). See the SVD detailed description for a description of condition variables.

◆ backsub() [2/2]

template<int Size = Dynamic, typename Precision = double>
template<int R, int C, typename P , typename B >
Matrix<Size,C, Precision> TooN::SymEigen< Size, Precision >::backsub ( const Matrix< R, C, P, B > &  rhs) const
inline

Calculate result of multiplying the (pseudo-)inverse of M by another matrix.

For a matrix \(A\), this calculates \(M^{\dagger}A\) by back substitution (i.e. without explictly calculating the (pseudo-)inverse). See the SVD detailed description for a description of condition variables.

◆ get_evalues()

template<int Size = Dynamic, typename Precision = double>
Vector<Size, Precision>& TooN::SymEigen< Size, Precision >::get_evalues ( )
inline

Returns the eigenvalues of the matrix.

The eigenvalues are listed in order, from the smallest to the largest. These are also the diagonal values of the matrix \(\Lambda\).

◆ get_evectors()

template<int Size = Dynamic, typename Precision = double>
Matrix<Size,Size,Precision>& TooN::SymEigen< Size, Precision >::get_evectors ( )
inline

Returns the eigenvectors of the matrix.

This returns \(U^T\), so that the rows of the matrix are the eigenvectors, which can be extracted using usual Matrix::operator[]() subscript operator. They are returned in order of the size of the corresponding eigenvalue, i.e. the vector with the largest eigenvalue is first.

◆ get_inv_diag()

template<int Size = Dynamic, typename Precision = double>
Vector<Size, Precision> TooN::SymEigen< Size, Precision >::get_inv_diag ( const double  condition) const
inline

Calculates the reciprocals of the eigenvalues of the matrix.

The vector invdiag lists the eigenvalues in order, from the largest (i.e. smallest reciprocal) to the smallest. These are also the diagonal values of the matrix \(Lambda^{-1}\). Any eigenvalues which are too small are set to zero (see the SVD detailed description for a description of the and condition variables).

◆ get_isqrtm()

template<int Size = Dynamic, typename Precision = double>
Matrix<Size, Size, Precision> TooN::SymEigen< Size, Precision >::get_isqrtm ( const double  condition = Internal::symeigen_condition_no) const
inline

Calculate the inverse square root of a matrix which is a matrix M such that M.T*M=A^-1.

Any square-rooted eigenvalues which are too small are set to zero (see the SVD detailed description for a description of the condition variables).

◆ get_pinv()

template<int Size = Dynamic, typename Precision = double>
Matrix<Size, Size, Precision> TooN::SymEigen< Size, Precision >::get_pinv ( const double  condition = Internal::symeigen_condition_no) const
inline

Calculate (pseudo-)inverse of the matrix.

This is not usually needed: if you need the inverse just to multiply it by a matrix or a vector, use one of the backsub() functions, which will be faster. See the SVD detailed description for a description of the pseudo-inverse and condition variables.

◆ get_sqrtm()

template<int Size = Dynamic, typename Precision = double>
Matrix<Size, Size, Precision> TooN::SymEigen< Size, Precision >::get_sqrtm ( ) const
inline

Calculate the square root of a matrix which is a matrix M such that M.T*M=A.


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