Low-level interface

The user will find here the low-level linear algebra interface. Most of the functions here are wrappers around the CBLAS or LAPACK API and the user should be careful when using them, see Overload of some BLAS and LAPACK functions to understand why.

lpk2mat

template<typename T, typename S>
matrix<T> castor::lpk2mat(std::vector<S> const &V, std::size_t m, std::size_t n)

See mat2lpk.

mat2lpk

template<typename T, typename S>
std::vector<T> castor::mat2lpk(matrix<S> const &A, std::size_t L)

See lpk2mat.

tgeev

template<typename T>
void castor::tgeev(std::string typ, int &n, std::vector<T> &A, std::vector<T> &E, std::vector<T> &V)
void castor::xgeev(char &jobl, char &jobr, int &n, std::vector<zlpk> &A, std::vector<zlpk> &E, std::vector<zlpk> &V, zlpk &wkopt, int &lwork, int &info)
void castor::xgeev(char &jobl, char &jobr, int &n, std::vector<clpk> &A, std::vector<clpk> &E, std::vector<clpk> &V, clpk &wkopt, int &lwork, int &info)

Eigenvalues and the left and/or right eigenvectors.

TGEEV computes for an N-by-N nonsymmetric matrix A of type T, the eigenvalues and, optionally, the left and/or right eigenvectors.

The left eigenvectors of A are the same as the right eigenvectors of A^t. If u(j) and v(j) are the left and right eigenvectors, respectively, corresponding to the eigenvalue lambda(j): (u(j)^t)*A = lambda(j)*(u(j)^t), A*v(j) = lambda(j) * v(j).

The computed eigenvectors are normalized to have Euclidean norm equal to 1 and largest component real.

tgels

template<typename T>
void castor::tgels(int &m, int &n, int &nrhs, std::vector<T> &A, std::vector<T> &B)
void castor::xgels(char &t, int &m, int &n, int &nrhs, std::vector<zlpk> &A, std::vector<zlpk> &B, int &l, zlpk &wkopt, int &lwork, int &info)
void castor::xgels(char &t, int &m, int &n, int &nrhs, std::vector<double> &A, std::vector<double> &B, int &l, double &wkopt, int &lwork, int &info)
void castor::xgels(char &t, int &m, int &n, int &nrhs, std::vector<clpk> &A, std::vector<clpk> &B, int &l, clpk &wkopt, int &lwork, int &info)
void castor::xgels(char &t, int &m, int &n, int &nrhs, std::vector<float> &A, std::vector<float> &B, int &l, float &wkopt, int &lwork, int &info)

Solve overdetermined or underdetermined linear systems.

TGELS solves overdetermined or underdetermined linear systems involving an M-by-N matrix A of type T, or its transpose, using a QR or LQ factorization of A. It is assumed that A has full rank.

If m >= n, it find the least squares solution of an overdetermined system, solving the least squares problem: minimize || B - A*X ||.

If m < n, it find the minimum norm solution of an underdetermined system A*X = B.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

tgemm

template<typename R, typename S>
void castor::tgemm(R alpha, matrix<double> const &A, matrix<double> const &B, S beta, matrix<double> &C)
template<typename R, typename S>
void castor::tgemm(R alpha, matrix<std::complex<float>> const &A, matrix<std::complex<float>> const &B, S beta, matrix<std::complex<float>> &C)
template<typename R, typename S>
void castor::tgemm(R alpha, matrix<std::complex<double>> const &A, matrix<std::complex<double>> const &B, S beta, matrix<std::complex<double>> &C)
template<typename R, typename S>
void castor::tgemm(R alpha, matrix<float> const &A, matrix<float> const &B, S beta, matrix<float> &C)

Matrix-matrix operations C = alpha*A*B + beta*C.

TGEMM performs one of the matrix-matrix operations using BLAS interface C = alpha*A*B + beta*C, alpha and beta are scalars of any type, and A, B and C are matrices of same type X, with A an m by k matrix, B a k by n matrix and C an m by n matrix.

See tgemm.

tgeqrf

template<typename T>
void castor::tgeqrf(int &m, int &n, std::vector<T> &A, std::vector<T> &R)
void castor::xgeqrf(int &m, int &n, int &l, std::vector<zlpk> &A, std::vector<zlpk> &tau, zlpk &wkopt, int &lwork, int &info)
void castor::xgeqrf(int &m, int &n, int &l, std::vector<double> &A, std::vector<double> &tau, double &wkopt, int &lwork, int &info)
void castor::xgeqrf(int &m, int &n, int &l, std::vector<clpk> &A, std::vector<clpk> &tau, clpk &wkopt, int &lwork, int &info)
void castor::xgeqrf(int &m, int &n, int &l, std::vector<float> &A, std::vector<float> &tau, float &wkopt, int &lwork, int &info)

Computes a QR factorization of a M-by-N matrix A of type T (A=Q*R).

tgesdd

template<typename T, typename T2>
void castor::tgesdd(std::string typ, int &m, int &n, std::vector<T> &A, std::vector<T2> &S, std::vector<T> &U, std::vector<T> &V)
void castor::xgesdd(char &job, int &m, int &n, int &l, std::vector<zlpk> &A, std::vector<double> &S, std::vector<zlpk> &U, std::vector<zlpk> &V, zlpk &wkopt, int &lwork, std::vector<int> &iwork, int &info)
void castor::xgesdd(char &job, int &m, int &n, int &l, std::vector<double> &A, std::vector<double> &S, std::vector<double> &U, std::vector<double> &V, double &wkopt, int &lwork, std::vector<int> &iwork, int &info)
void castor::xgesdd(char &job, int &m, int &n, int &l, std::vector<clpk> &A, std::vector<float> &S, std::vector<clpk> &U, std::vector<clpk> &V, clpk &wkopt, int &lwork, std::vector<int> &iwork, int &info)
void castor::xgesdd(char &job, int &m, int &n, int &l, std::vector<float> &A, std::vector<float> &S, std::vector<float> &U, std::vector<float> &V, float &wkopt, int &lwork, std::vector<int> &iwork, int &info)

Singular value decomposition of a rectangular matrix A.

TGESDD computes the singular value decomposition (SVD) of a m-by-n matrix A of type T, optionally computing the left and/or right singular vectors. If singular vectors are desired, it uses a divide and conquer algorithm.

The SVD is written as : A = U*SIGMA*VT where SIGMA is an m-by-n matrix which is zero except for its min(m,n) diagonal elements, U is an m-by-m orthogonal matrix and VT (V transposed) is an n-by-n orthogonal matrix. The diagonal elements of SIGMA are the singular values of A; they are real and non-negative, and are returned in descending order. The first min(m, n) columns of U and V are the left and right singular vectors of A.

Note that the routine returns VT, not V.

tgesv

template<typename T>
void castor::tgesv(int &n, int &nrhs, std::vector<T> &A, std::vector<T> &B)
void castor::xgesv(int &n, int &nrhs, std::vector<zlpk> &A, std::vector<zlpk> &B, std::vector<int> &ipiv, int &info)
void castor::xgesv(int &n, int &nrhs, std::vector<double> &A, std::vector<double> &B, std::vector<int> &ipiv, int &info)
void castor::xgesv(int &n, int &nrhs, std::vector<clpk> &A, std::vector<clpk> &B, std::vector<int> &ipiv, int &info)
void castor::xgesv(int &n, int &nrhs, std::vector<float> &A, std::vector<float> &B, std::vector<int> &ipiv, int &info)

Solve linear system equations A*X=B.

TGESV computes the solution to a system of linear equations A*X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices, both of type T.

The LU decomposition with partial pivoting and row interchanges is used to factor A as: A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

tgetrf

template<typename T, typename S>
void castor::tgetrf(int &m, int &n, std::vector<T> &A, std::vector<T> &P, matrix<S> &U)
void castor::xgetrf(int &m, int &n, std::vector<zlpk> &A, std::vector<int> &ipiv, int &info)
void castor::xgetrf(int &m, int &n, std::vector<double> &A, std::vector<int> &ipiv, int &info)
void castor::xgetrf(int &m, int &n, std::vector<clpk> &A, std::vector<int> &ipiv, int &info)
void castor::xgetrf(int &m, int &n, std::vector<float> &A, std::vector<int> &ipiv, int &info)

LU factorization using partial pivoting.

TGETRF computes an LU factorization of a general M-by-N matrix using partial pivoting with row interchanges. The factorization has the form: P * A = L * U where P is a row-pivot matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.