FILTLAN  1.0a
Functions
Eigenvalues / Eigenvectors of a Symmetric Tridiaognal Matrix

Functions

int SymmetricTridiagoanlEigenSolver (Vector &eigVal, mkIndex n, const Real *diag, const Real *sdiag)
 Compute all eigenvalues (but not eigenvectors) of a symmetric tridiagonal matrix T. This function invokes the LAPACK routine dstev_() (if Real is double) or stev_() (if Real is float). More...
 
int SymmetricTridiagoanlEigenSolver (Vector &eigVal, Matrix &eigVec, mkIndex n, const Real *diag, const Real *sdiag)
 Compute all eigenvalues and eigenvectors of a symmetric tridiagonal matrix T. This function invokes the LAPACK routine dstev_() (if Real is double) or stev_() (if Real is float). More...
 
void reportTroubleIfAny (std::ostream &outerr, int info, mkIndex n=0)
 Print out an error message to std::ostream outerr if info $\neq$ 0 which signifies an error that occurred in the LAPACK routine xSTEV. More...
 

Detailed Description

Function Documentation

◆ reportTroubleIfAny()

void reportTroubleIfAny ( std::ostream &  outerr,
int  info,
mkIndex  n = 0 
)

Print out an error message to std::ostream outerr if info $\neq$ 0 which signifies an error that occurred in the LAPACK routine xSTEV.

Parameters
outerris an std::ostream for printing out the error message when info $\neq$ 0.
infois the flag returned from the LAPACK routine dstev_() (if Real is double) or sstev_() (if Real is float).
  • =0: successful exit (nothing will be printed).
  • <0: if info=-i, the ith argument had an illegal value.
  • >0: if info=i, the algorithm failed to converge; i off-diagonal elements did not converge to zero.
nis the dimension of the symmetric tridiagonal matrix. This information will also be printed out if n is provided (i.e. n $\neq$ 0).
See also
SymmetricTridiagoanlEigenSolver().

◆ SymmetricTridiagoanlEigenSolver() [1/2]

int SymmetricTridiagoanlEigenSolver ( Vector &  eigVal,
mkIndex  n,
const Real *  diag,
const Real *  sdiag 
)

Compute all eigenvalues (but not eigenvectors) of a symmetric tridiagonal matrix T. This function invokes the LAPACK routine dstev_() (if Real is double) or stev_() (if Real is float).

Parameters
nis the dimension of the symmatrix tridiagonal matrix T.
diag[],sdiag[]define the symmetric tridiagonal matrix T, whose the diagonal elements are diag[0,...,n-1] and subdiagonal elements are sdiag[0,...,n-2]. In other words,
  • T(i,i) = diag[i-1] for i = 1,...,n.
  • T(i+1,i) = T(i,i+1) = sdiag[i-1] for i = 1,...n-1.
eigValis the output vector of length n containing all eigenvalues of T in ascending order.
Returns
The flag returned by the LAPACK routine dstev_() (if Real is double) or stev_() (if Real is float).
See also
reportTroubleIfAny().

◆ SymmetricTridiagoanlEigenSolver() [2/2]

int SymmetricTridiagoanlEigenSolver ( Vector &  eigVal,
Matrix &  eigVec,
mkIndex  n,
const Real *  diag,
const Real *  sdiag 
)

Compute all eigenvalues and eigenvectors of a symmetric tridiagonal matrix T. This function invokes the LAPACK routine dstev_() (if Real is double) or stev_() (if Real is float).

Parameters
nis the dimension of the symmatrix tridiagonal matrix T.
diag[],sdiag[]define the symmetric tridiagonal matrix; the diagonal elements are diag[0,...,n-1] and subdiagonal elements are sdiag[0,...,n-2]. In other words,
  • T(i,i) = diag[i-1] for i = 1,...,n.
  • T(i+1,i) = T(i,i+1) = sdiag[i-1] for i = 1,...n-1.
eigValis the output vector of length n containing all eigenvalues of T in ascending order.
eigVecis the output n-by-n matrix with columns as eigenvectors of T, sorted with respect to eigVal.
Returns
The flag returned by dstev_() (if Real is double) or stev_() (if Real is float).
See also
reportTroubleIfAny().