Namespace xlifepp::internalEigenSolver#

namespace internalEigenSolver

namespace dedicated to internal eigensolver of XLiFE++

Functions

template<typename MatrixType, typename VectorsType, typename CoeffsType>
void applyBlockHouseholderOnTheLeft(MatrixType &mat, VectorsType &vectors, const CoeffsType &hCoeffs)
template<typename MatrixType>
MatrixType directSwapping(MatrixType &matT, unsigned int p, unsigned int q, real_t gamma)
template<typename MatrixType>
void doComputeEigenvectorsComplexSolverInPlace(const real_t matrixnorm, const MatrixType &schurT, const MatrixType &schurU, MatrixType &eigVec)
template<typename MatrixType>
void doComputeEigenvectorsRealSolverInPlace(const real_t matrixnorm, const MatrixType &schurT, const MatrixType &schurU, MatrixType &eigVec)
template<typename MatrixQR, typename HCoeffs>
void householderQRinplaceBlocked(MatrixQR &mat, HCoeffs &hCoeffs, Index maxBlockSize = 32)
template<typename MatrixQR, typename HCoeffs>
void householderQRinplaceUnblocked(MatrixQR &mat, HCoeffs &hCoeffs)
template<typename TriangularFactorType, typename VectorsType, typename CoeffsType>
void makeBlockHouseholderTriangularFactor(TriangularFactorType &triFactor, VectorsType &vectors, const CoeffsType &hCoeffs)
inline void printOutDebugInfoEigenProblem(const string_t &nameOfClass, const string_t &message)
template<typename MatrixType>
void swapComplexSchurInPlace(MatrixType &matrixT, const std::vector<int> &order)
template<typename MatrixType>
void swapComplexSchurInPlace(MatrixType &matrixT, MatrixType &matrixQ, const std::vector<int> &order)
template<typename MatrixType>
void swapRealSchurInPlace(MatrixType &matrixT, MatrixType &matrixQ, unsigned int ifst, unsigned int ilst)

Given the quasi-triangular matrix T and orthogonal matrix Q obtained from the real Schur decomposition, this function reorders the eigenvalues appearing on the (block) diagonal of matrix T The diagonal element of block of T with row index ifst is moved to row ilst.

This is implemented from “On Swapping Diagonal Blocks in Real Schur Form”

Parameters:
  • matrixT[inout] quasi-triangular matrix of real Schur decomposition

  • matrixQ[inout] orthogonal of real Schur decomposition

  • ifst[in] index of row needs moving

  • ilst[in] index of destination row

inline void testErrorEigenProblem(bool condition, const string_t &s)
inline void testErrorEigenProblemMultVec(bool condition, const string_t &s)
inline void testWarningEigenProblem(bool condition, const string_t &s)
template<typename MatrixType, typename DiagonalType, typename SubDiagonalType>
void tridiagonalizationInplace(MatrixType &mat, DiagonalType &diag, SubDiagonalType &subdiag, bool extractQ)

Performs a full tridiagonalization in place.

Computes the tridiagonal decomposition of the selfadjoint matrix mat in place such that \( mat = Q T Q^* \) where \( Q \) is unitary and \( T \) a real symmetric tridiagonal matrix.

The tridiagonal matrix T is passed to the output parameters diag and subdiag. If extractQ is true, then the orthogonal matrix Q is passed to mat. Otherwise the lower part of the matrix mat is destroyed.

The vectors diag and subdiag are not resized. The function assumes that they are already of the correct size. The length of the vector diag should equal the number of rows in mat, and the length of the vector subdiag should be one left.

See also

class Tridiagonalization

Parameters:
  • mat[inout] On input, the selfadjoint matrix whose tridiagonal decomposition is to be computed. Only the lower triangular part referenced. The rest is left unchanged. On output, the orthogonal matrix Q in the decomposition if extractQ is true.

  • diag[out] The diagonal of the tridiagonal matrix T in the decomposition.

  • subdiag[out] The subdiagonal of the tridiagonal matrix T in the decomposition.

  • extractQ[in] If true, the orthogonal matrix Q in the decomposition is computed and stored in mat.

Note

Currently, it requires two temporary vectors to hold the intermediate Householder coefficients, and to reconstruct the matrix Q from the Householder reflectors.

template<typename MatrixType, typename CoeffVectorType>
void tridiagonalizationInplace(MatrixType &matA, CoeffVectorType &hCoeffs)

Performs a tridiagonal decomposition of the selfadjoint matrix matA in-place.

On output, the tridiagonal selfadjoint matrix T is stored in the diagonal and lower sub-diagonal of the matrix matA. The unitary matrix Q is represented in a compact way as a product of Householder reflectors \( H_i \) such that: \( Q = H_{N-1} \ldots H_1 H_0 \). The Householder reflectors are defined as \( H_i = (I - h_i v_i v_i^T) \) where \( h_i = hCoeffs[i]\) is the \( i \)th Householder coefficient and \( v_i \) is the Householder vector defined by \( v_i = [ 0, \ldots, 0, 1, matA(i+2,i), \ldots, matA(N-1,i) ]^T \).

Implemented from Golub’s “Matrix Computations”, algorithm 8.3.1.

Parameters:
  • matA[inout] On input the selfadjoint matrix. Only the lower triangular part is referenced. On output, the strict upper part is left unchanged, and the lower triangular part represents the T and Q matrices in packed format has detailed below.

  • hCoeffs[out] returned Householder coefficients (see below)

template<typename MatrixType, typename CoeffVectorType>
void tridiagonalizationUpperInplace(MatrixType &matA, CoeffVectorType &hCoeffs)
template<typename MatrixType, bool IsComplex>
struct ComplexSchurReduceToHessenberg
template<typename MatrixType, bool IsComplex>
struct ComputeSchurForm
struct EigenSolverVerbose
#include <EigenTestError.hpp>

struct to handle verbose level

template<typename MatrixType, bool IsComplex>
struct EigenValueComputation
template<typename MatrixType, bool IsComplex>
struct EigenVectorComputation
template<typename MatrixType, bool IsComplex>
struct SortSchurForm