LAPACK Functions

Table of Contents

The Blaze library makes extensive use of the LAPACK functionality for various compute tasks (including the decomposition, inversion and the computation of the determinant of dense matrices). For this purpose, Blaze implements several convenient C++ wrapper functions for all required LAPACK functions. The following sections give a complete overview of all available LAPACK wrapper functions. For more details on the individual LAPACK functions see the Blaze function documentation or the LAPACK online documentation browser:

http://www.netlib.org/lapack/explore-html/

Note
All functions only work for general, non-adapted matrices with float, double, complex<float>, or complex<double> element type. The attempt to call the function with adaptors or matrices of any other element type results in a compile time error!
All functions can only be used if the fitting LAPACK library is available and linked to the final executable. Otherwise a call to this function will result in a linker error.
For performance reasons all functions do only provide the basic exception safety guarantee, i.e. in case an exception is thrown the given matrix may already have been modified.


Matrix Decomposition


The following functions decompose/factorize the given dense matrix. Based on this decomposition the matrix can be inverted or used to solve a linear system of equations.


LU Decomposition

The following functions provide an interface for the LAPACK functions sgetrf(), dgetrf(), cgetrf(), and zgetrf(), which compute the LU decomposition for the given general matrix:

namespace blaze {
void getrf( int m, int n, float* A, int lda, int* ipiv, int* info );
void getrf( int m, int n, double* A, int lda, int* ipiv, int* info );
void getrf( int m, int n, complex<float>* A, int lda, int* ipiv, int* info );
void getrf( int m, int n, complex<double>* A, int lda, int* ipiv, int* info );
template< typename MT, bool SO >
void getrf( DenseMatrix<MT,SO>& A, int* ipiv );
} // namespace blaze

The decomposition has the form

\[ A = P \cdot L \cdot U, \]


where P is a permutation matrix, L is a lower unitriangular matrix, and U is an upper triangular matrix. The resulting decomposition is stored within A: In case of a column-major matrix, L is stored in the lower part of A and U is stored in the upper part. The unit diagonal elements of L are not stored. In case A is a row-major matrix the result is transposed.

Note
The LU decomposition will never fail, even for singular matrices. However, in case of a singular matrix the resulting decomposition cannot be used for a matrix inversion or solving a linear system of equations.


LDLT Decomposition

The following functions provide an interface for the LAPACK functions ssytrf(), dsytrf(), csytrf(), and zsytrf(), which compute the LDLT (Bunch-Kaufman) decomposition for the given symmetric indefinite matrix:

namespace blaze {
void sytrf( char uplo, int n, float* A, int lda, int* ipiv, float* work, int lwork, int* info );
void sytrf( char uplo, int n, double* A, int lda, int* ipiv, double* work, int lwork, int* info );
void sytrf( char uplo, int n, complex<float>* A, int lda, int* ipiv, complex<float>* work, int lwork, int* info );
void sytrf( char uplo, int n, complex<double>* A, int lda, int* ipiv, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void sytrf( DenseMatrix<MT,SO>& A, char uplo, int* ipiv );
} // namespace blaze

The decomposition has the form

\[ A = U D U^{T} \texttt{ (if uplo = 'U'), or } A = L D L^{T} \texttt{ (if uplo = 'L'), } \]

where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L' the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U' the result is stored in the upper part and the lower part remains untouched.

Note
The Bunch-Kaufman decomposition will never fail, even for singular matrices. However, in case of a singular matrix the resulting decomposition cannot be used for a matrix inversion or solving a linear system of equations.


LDLH Decomposition

The following functions provide an interface for the LAPACK functions chetrf() and zsytrf(), which compute the LDLH (Bunch-Kaufman) decomposition for the given Hermitian indefinite matrix:

namespace blaze {
void hetrf( char uplo, int n, complex<float>* A, int lda, int* ipiv, complex<float>* work, int lwork, int* info );
void hetrf( char uplo, int n, complex<double>* A, int lda, int* ipiv, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void hetrf( DenseMatrix<MT,SO>& A, char uplo, int* ipiv );
} // namespace blaze

The decomposition has the form

\[ A = U D U^{H} \texttt{ (if uplo = 'U'), or } A = L D L^{H} \texttt{ (if uplo = 'L'), } \]

where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L' the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U' the result is stored in the upper part and the lower part remains untouched.

Note
The Bunch-Kaufman decomposition will never fail, even for singular matrices. However, in case of a singular matrix the resulting decomposition cannot be used for a matrix inversion or solving a linear system of equations.


Cholesky Decomposition

The following functions provide an interface for the LAPACK functions spotrf(), dpotrf(), cpotrf(), and zpotrf(), which compute the Cholesky (LLH) decomposition for the given positive definite matrix:

namespace blaze {
void potrf( char uplo, int n, float* A, int lda, int* info );
void potrf( char uplo, int n, double* A, int lda, int* info );
void potrf( char uplo, int n, complex<float>* A, int lda, int* info );
void potrf( char uplo, int n, complex<double>* A, int lda, int* info );
template< typename MT, bool SO >
void potrf( DenseMatrix<MT,SO>& A, char uplo );
} // namespace blaze

The decomposition has the form

\[ A = U^{T} U \texttt{ (if uplo = 'U'), or } A = L L^{T} \texttt{ (if uplo = 'L'), } \]

where U is an upper triangular matrix and L is a lower triangular matrix. The Cholesky decomposition fails if the given matrix A is not a positive definite matrix. In this case a std::std::invalid_argument exception is thrown.


QR Decomposition

The following functions provide an interface for the LAPACK functions sgeqrf(), dgeqrf(), cgeqrf(), and zgeqrf(), which compute the QR decomposition of the given general matrix:

namespace blaze {
void geqrf( int m, int n, float* A, int lda, float* tau, float* work, int lwork, int* info );
void geqrf( int m, int n, double* A, int lda, double* tau, double* work, int lwork, int* info );
void geqrf( int m, int n, complex<float>* A, int lda, complex<float>* tau, complex<float>* work, int lwork, int* info );
void geqrf( int m, int n, complex<double>* A, int lda, complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void geqrf( DenseMatrix<MT,SO>& A, typename MT::ElementType* tau );
} // namespace blaze

The decomposition has the form

\[ A = Q \cdot R, \]

where the Q is represented as a product of elementary reflectors

\[ Q = H(1) H(2) . . . H(k) \texttt{, with k = min(m,n).} \]

Each H(i) has the form

\[ H(i) = I - tau \cdot v \cdot v^T, \]

where tau is a real scalar, and v is a real vector with v(0:i-1) = 0 and v(i) = 1. v(i+1:m) is stored on exit in A(i+1:m,i), and tau in tau(i). Thus on exit the elements on and above the diagonal of the matrix contain the min(m,n)-by-n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array tau, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors.

The following functions provide an interface for the LAPACK functions sorgqr(), dorgqr(), cungqr(), and zunqqr(), which reconstruct the Q matrix from a QR decomposition:

namespace blaze {
void orgqr( int m, int n, int k, float* A, int lda, const float* tau, float* work, int lwork, int* info );
void orgqr( int m, int n, int k, double* A, int lda, const double* tau, double* work, int lwork, int* info );
void ungqr( int m, int n, int k, complex<float>* A, int lda, const complex<float>* tau, complex<float>* work, int lwork, int* info );
void ungqr( int m, int n, int k, complex<double>* A, int lda, const complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void orgqr( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
template< typename MT, bool SO >
void ungqr( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
} // namespace blaze


RQ Decomposition

The following functions provide an interface for the LAPACK functions sgerqf(), dgerqf(), cgerqf(), and zgerqf(), which compute the RQ decomposition of the given general matrix:

namespace blaze {
void gerqf( int m, int n, float* A, int lda, float* tau, float* work, int lwork, int* info );
void gerqf( int m, int n, double* A, int lda, double* tau, double* work, int lwork, int* info );
void gerqf( int m, int n, complex<float>* A, int lda, complex<float>* tau, complex<float>* work, int lwork, int* info );
void gerqf( int m, int n, complex<double>* A, int lda, complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void gerqf( DenseMatrix<MT,SO>& A, typename MT::ElementType* tau );
} // namespace blaze

The decomposition has the form

\[ A = R \cdot Q, \]

where the Q is represented as a product of elementary reflectors

\[ Q = H(1) H(2) . . . H(k) \texttt{, with k = min(m,n).} \]

Each H(i) has the form

\[ H(i) = I - tau \cdot v \cdot v^T, \]

where tau is a real scalar, and v is a real vector with v(n-k+i+1:n) = 0 and v(n-k+i) = 1. v(1:n-k+i-1) is stored on exit in A(m-k+i,1:n-k+i-1), and tau in tau(i). Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R; the remaining elements in combination with the array tau represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors.

The following functions provide an interface for the LAPACK functions sorgrq(), dorgrq(), cungrq(), and zunqrq(), which reconstruct the Q matrix from a RQ decomposition:

namespace blaze {
void orgrq( int m, int n, int k, float* A, int lda, const float* tau, float* work, int lwork, int* info );
void orgrq( int m, int n, int k, double* A, int lda, const double* tau, double* work, int lwork, int* info );
void ungrq( int m, int n, int k, complex<float>* A, int lda, const complex<float>* tau, complex<float>* work, int lwork, int* info );
void ungrq( int m, int n, int k, complex<double>* A, int lda, const complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void orgrq( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
template< typename MT, bool SO >
void ungrq( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
} // namespace blaze


QL Decomposition

The following functions provide an interface for the LAPACK functions sgeqlf(), dgeqlf(), cgeqlf(), and zgeqlf(), which compute the QL decomposition of the given general matrix:

namespace blaze {
void geqlf( int m, int n, float* A, int lda, float* tau, float* work, int lwork, int* info );
void geqlf( int m, int n, double* A, int lda, double* tau, double* work, int lwork, int* info );
void geqlf( int m, int n, complex<float>* A, int lda, complex<float>* tau, complex<float>* work, int lwork, int* info );
void geqlf( int m, int n, complex<double>* A, int lda, complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void geqlf( DenseMatrix<MT,SO>& A, typename MT::ElementType* tau );
} // namespace blaze

The decomposition has the form

\[ A = Q \cdot L, \]

where the Q is represented as a product of elementary reflectors

\[ Q = H(k) . . . H(2) H(1) \texttt{, with k = min(m,n).} \]

Each H(i) has the form

\[ H(i) = I - tau \cdot v \cdot v^T, \]

where tau is a real scalar, and v is a real vector with v(m-k+i+1:m) = 0 and v(m-k+i) = 1. v(1:m-k+i-1) is stored on exit in A(1:m-k+i-1,n-k+i), and tau in tau(i). Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L; the remaining elements in combination with the array tau represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors.

The following functions provide an interface for the LAPACK functions sorgql(), dorgql(), cungql(), and zunqql(), which reconstruct the Q matrix from an QL decomposition:

namespace blaze {
void orgql( int m, int n, int k, float* A, int lda, const float* tau, float* work, int lwork, int* info );
void orgql( int m, int n, int k, double* A, int lda, const double* tau, double* work, int lwork, int* info );
void ungql( int m, int n, int k, complex<float>* A, int lda, const complex<float>* tau, complex<float>* work, int lwork, int* info );
void ungql( int m, int n, int k, complex<double>* A, int lda, const complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void orgql( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
template< typename MT, bool SO >
void ungql( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
} // namespace blaze


LQ Decomposition

The following functions provide an interface for the LAPACK functions sgelqf(), dgelqf(), cgelqf(), and zgelqf(), which compute the LQ decomposition of the given general matrix:

namespace blaze {
void gelqf( int m, int n, float* A, int lda, float* tau, float* work, int lwork, int* info );
void gelqf( int m, int n, double* A, int lda, double* tau, double* work, int lwork, int* info );
void gelqf( int m, int n, complex<float>* A, int lda, complex<float>* tau, complex<float>* work, int lwork, int* info );
void gelqf( int m, int n, complex<double>* A, int lda, complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void gelqf( DenseMatrix<MT,SO>& A, typename MT::ElementType* tau );
} // namespace blaze

The decomposition has the form

\[ A = L \cdot Q, \]

where the Q is represented as a product of elementary reflectors

\[ Q = H(k) . . . H(2) H(1) \texttt{, with k = min(m,n).} \]

Each H(i) has the form

\[ H(i) = I - tau \cdot v \cdot v^T, \]

where tau is a real scalar, and v is a real vector with v(0:i-1) = 0 and v(i) = 1. v(i+1:n) is stored on exit in A(i,i+1:n), and tau in tau(i). Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tau, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors.

The following functions provide an interface for the LAPACK functions sorglq(), dorglq(), cunglq(), and zunqlq(), which reconstruct the Q matrix from an LQ decomposition:

namespace blaze {
void orglq( int m, int n, int k, float* A, int lda, const float* tau, float* work, int lwork, int* info );
void orglq( int m, int n, int k, double* A, int lda, const double* tau, double* work, int lwork, int* info );
void unglq( int m, int n, int k, complex<float>* A, int lda, const complex<float>* tau, complex<float>* work, int lwork, int* info );
void unglq( int m, int n, int k, complex<double>* A, int lda, const complex<double>* tau, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void orglq( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
template< typename MT, bool SO >
void unglq( DenseMatrix<MT,SO>& A, const typename MT::ElementType* tau );
} // namespace blaze


Matrix Inversion


Given a matrix that has already been decomposed, the following functions can be used to invert the matrix in-place.


LU-based Inversion

The following functions provide an interface for the LAPACK functions sgetri(), dgetri(), cgetri(), and zgetri(), which invert a general matrix that has already been decomposed by an LU Decomposition :

namespace blaze {
void getri( int n, float* A, int lda, const int* ipiv, float* work, int lwork, int* info );
void getri( int n, double* A, int lda, const int* ipiv, double* work, int lwork, int* info );
void getri( int n, complex<float>* A, int lda, const int* ipiv, complex<float>* work, int lwork, int* info );
void getri( int n, complex<double>* A, int lda, const int* ipiv, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO >
void getri( DenseMatrix<MT,SO>& A, const int* ipiv );
} // namespace blaze

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


LDLT-based Inversion

The following functions provide an interface for the LAPACK functions ssytri(), dsytri(), csytri(), and zsytri(), which invert a symmetric indefinite matrix that has already been decomposed by an LDLT Decomposition :

namespace blaze {
void sytri( char uplo, int n, float* A, int lda, const int* ipiv, float* work, int* info );
void sytri( char uplo, int n, double* A, int lda, const int* ipiv, double* work, int* info );
void sytri( char uplo, int n, complex<float>* A, int lda, const int* ipiv, complex<float>* work, int* info );
void sytri( char uplo, int n, complex<double>* A, int lda, const int* ipiv, complex<double>* work, int* info );
template< typename MT, bool SO >
void sytri( DenseMatrix<MT,SO>& A, char uplo, const int* ipiv );
} // namespace blaze

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


LDLH-based Inversion

The following functions provide an interface for the LAPACK functions chetri() and zhetri(), which invert an Hermitian indefinite matrix that has already been decomposed by an LDLH Decomposition :

namespace blaze {
void hetri( char uplo, int n, complex<float>* A, int lda, const int* ipiv, complex<float>* work, int* info );
void hetri( char uplo, int n, complex<double>* A, int lda, const int* ipiv, complex<double>* work, int* info );
template< typename MT, bool SO >
void hetri( DenseMatrix<MT,SO>& A, char uplo, const int* ipiv );
} // namespace blaze

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


Cholesky-based Inversion

The following functions provide an interface for the LAPACK functions spotri(), dpotri(), cpotri(), and zpotri(), which invert a positive definite matrix that has already been decomposed by an Cholesky Decomposition :

namespace blaze {
void potri( char uplo, int n, float* A, int lda, int* info );
void potri( char uplo, int n, double* A, int lda, int* info );
void potri( char uplo, int n, complex<float>* A, int lda, int* info );
void potri( char uplo, int n, complex<double>* A, int lda, int* info );
template< typename MT, bool SO >
void potri( DenseMatrix<MT,SO>& A, char uplo );
} // namespace blaze

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


Inversion of Triangular Matrices

The following functions provide an interface for the LAPACK functions strtri(), dtrtri(), ctrtri(), and ztrtri(), which invert the given triangular matrix in-place:

namespace blaze {
void trtri( char uplo, char diag, int n, float* A, int lda, int* info );
void trtri( char uplo, char diag, int n, double* A, int lda, int* info );
void trtri( char uplo, char diag, int n, complex<float>* A, int lda, int* info );
void trtri( char uplo, char diag, int n, complex<double>* A, int lda, int* info );
template< typename MT, bool SO >
void trtri( DenseMatrix<MT,SO>& A, char uplo, char diag );
} // namespace blaze

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


Substitution


Given a matrix that has already been decomposed the following functions can be used to perform the forward/backward substitution step to compute the solution to a system of linear equations. Note that depending on the storage order of the system matrix and the given right-hand side the functions solve different equation systems:

Single right-hand side:

Multiple right-hand sides:

In this context the general system matrix A is a n-by-n matrix that has already been factorized by the according decomposition function, x and b are n-dimensional vectors and X and B are either row-major m-by-n matrices or column-major n-by-m matrices.


LU-based Substitution

The following functions provide an interface for the LAPACK functions sgetrs(), dgetrs(), cgetrs(), and zgetrs(), which perform the substitution step for a general matrix that has already been decomposed by an LU Decomposition :

namespace blaze {
void getrs( char trans, int n, int nrhs, const float* A, int lda, const int* ipiv, float* B, int ldb, int* info );
void getrs( char trans, int n, int nrhs, const double* A, int lda, const int* ipiv, double* B, int ldb, int* info );
void getrs( char trans, int n, const complex<float>* A, int lda, const int* ipiv, complex<float>* B, int ldb, int* info );
void getrs( char trans, int n, const complex<double>* A, int lda, const int* ipiv, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void getrs( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char trans, const int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void getrs( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char trans, const int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations. The function fails if ...

The first four functions report failure via the info argument, the last two functions throw a std::invalid_argument exception in case of an error.


LDLT-based Substitution

The following functions provide an interface for the LAPACK functions ssytrs(), dsytrs(), csytrs(), and zsytrs(), which perform the substitution step for a symmetric indefinite matrix that has already been decomposed by an LDLT Decomposition :

namespace blaze {
void sytrs( char uplo, int n, int nrhs, const float* A, int lda, const int* ipiv, float* B, int ldb, int* info );
void sytrs( char uplo, int n, int nrhs, const double* A, int lda, const int* ipiv, double* B, int ldb, int* info );
void sytrs( char uplo, int n, int nrhs, const complex<float>* A, int lda, const int* ipiv, complex<float>* B, int ldb, int* info );
void sytrs( char uplo, int n, int nrhs, const complex<double>* A, int lda, const int* ipiv, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void sytrs( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, const int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void sytrs( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo, const int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations. The function fails if ...

The first four functions report failure via the info argument, the last two functions throw a std::invalid_argument exception in case of an error.


LDLH-based Substitution

The following functions provide an interface for the LAPACK functions chetrs(), and zhetrs(), which perform the substitution step for an Hermitian indefinite matrix that has already been decomposed by an LDLH Decomposition :

namespace blaze {
void hetrs( char uplo, int n, int nrhs, const complex<float>* A, int lda, const int* ipiv, complex<float>* B, int ldb, int* info );
void hetrs( char uplo, int n, int nrhs, const complex<double>* A, int lda, const int* ipiv, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void hetrs( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, const int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void hetrs( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo, const int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations. The function fails if ...

The first two functions report failure via the info argument, the last two functions throw a std::invalid_argument exception in case of an error.


Cholesky-based Substitution

The following functions provide an interface for the LAPACK functions spotrs(), dpotrs(), cpotrs(), and zpotrs(), which perform the substitution step for a positive definite matrix that has already been decomposed by an Cholesky Decomposition :

namespace blaze {
void potrs( char uplo, int n, int nrhs, const float* A, int lda, float* B, int ldb, int* info );
void potrs( char uplo, int n, int nrhs, const double* A, int lda, double* B, int ldb, int* info );
void potrs( char uplo, int n, int nrhs, const complex<float>* A, int lda, complex<float>* B, int ldb, int* info );
void potrs( char uplo, int n, int nrhs, const complex<double>* A, int lda, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void potrs( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void potrs( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations. The function fails if ...

The first two functions report failure via the info argument, the last two functions throw a std::invalid_argument exception in case of an error.


Substitution for Triangular Matrices

The following functions provide an interface for the LAPACK functions strtrs(), dtrtrs(), ctrtrs(), and ztrtrs(), which perform the substitution step for a triangular matrix:

namespace blaze {
void trtrs( char uplo, char trans, char diag, int n, int nrhs, const float* A, int lda, float* B, int ldb, int* info );
void trtrs( char uplo, char trans, char diag, int n, int nrhs, const double* A, int lda, double* B, int ldb, int* info );
void trtrs( char uplo, char trans, char diag, int n, int nrhs, const complex<float>* A, int lda, complex<float>* B, int ldb, int* info );
void trtrs( char uplo, char trans, char diag, int n, int nrhs, const complex<double>* A, int lda, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void trtrs( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, char trans, char diag );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void trtrs( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo, char trans, char diag );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations. The function fails if ...

The first four functions report failure via the info argument, the last two functions throw a std::invalid_argument exception in case of an error.


Linear System Solver


The following functions represent compound functions that perform both the decomposition step as well as the substitution step to compute the solution to a system of linear equations. Note that depending on the storage order of the system matrix and the given right-hand side the functions solve different equation systems:

Single right-hand side:

Multiple right-hand sides:

In this context the general system matrix A is a n-by-n matrix that has already been factorized by the according decomposition function, x and b are n-dimensional vectors and X and B are either row-major m-by-n matrices or column-major n-by-m matrices.

LU-based Linear System Solver

The following functions provide an interface for the LAPACK functions sgesv(), dgesv(), cgesv(), and zgesv(), which combine an LU Decomposition and the according LU-based Substitution :

namespace blaze {
void gesv( int n, int nrhs, float* A, int lda, int* ipiv, float* B, int ldb, int* info );
void gesv( int n, int nrhs, double* A, int lda, int* ipiv, double* B, int ldb, int* info );
void gesv( int n, int nrhs, complex<float>* A, int lda, int* ipiv, complex<float>* B, int ldb, int* info );
void gesv( int n, int nrhs, complex<double>* A, int lda, int* ipiv, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void gesv( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void gesv( DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations and A has been decomposed by means of an LU Decomposition.

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


LDLT-based Linear System Solver

The following functions provide an interface for the LAPACK functions ssysv(), dsysv(), csysv(), and zsysv(), which combine an LDLT Decomposition and the according LDLT-based Substitution :

namespace blaze {
void sysv( char uplo, int n, int nrhs, float* A, int lda, int* ipiv, float* B, int ldb, float* work, int lwork, int* info );
void sysv( char uplo, int n, int nrhs, double* A, int lda, int* ipiv, double* B, int ldb, double* work, int lwork, int* info );
void sysv( char uplo, int n, int nrhs, complex<float>* A, int lda, int* ipiv, complex<float>* B, int ldb, complex<float>* work, int lwork, int* info );
void sysv( char uplo, int n, int nrhs, complex<double>* A, int lda, int* ipiv, complex<double>* B, int ldb, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void sysv( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void sysv( DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo, int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations and A has been decomposed by means of an LDLT Decomposition.

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


LDLH-based Linear System Solver

The following functions provide an interface for the LAPACK functions shesv(), dhesv(), chesv(), and zhesv(), which combine an LDLH Decomposition and the according LDLH-based Substitution :

namespace blaze {
void hesv( char uplo, int n, int nrhs, complex<float>* A, int lda, int* ipiv, complex<float>* B, int ldb, complex<float>* work, int lwork, int* info );
void hesv( char uplo, int n, int nrhs, complex<double>* A, int lda, int* ipiv, complex<double>* B, int ldb, complex<double>* work, int lwork, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void hesv( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, int* ipiv );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void hesv( DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo, int* ipiv );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations and A has been decomposed by means of an LDLH Decomposition.

The functions fail if ...

The first two functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


Cholesky-based Linear System Solver

The following functions provide an interface for the LAPACK functions sposv(), dposv(), cposv(), and zposv(), which combine an Cholesky Decomposition and the according Cholesky-based Substitution :

namespace blaze {
void posv( char uplo, int n, int nrhs, float* A, int lda, float* B, int ldb, int* info );
void posv( char uplo, int n, int nrhs, double* A, int lda, double* B, int ldb, int* info );
void posv( char uplo, int n, int nrhs, complex<float>* A, int lda, complex<float>* B, int ldb, int* info );
void posv( char uplo, int n, int nrhs, complex<double>* A, int lda, complex<double>* B, int ldb, int* info );
template< typename MT, bool SO, typename VT, bool TF >
void posv( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo );
template< typename MT1, bool SO1, typename MT2, bool SO2 >
void posv( DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B, char uplo );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations and A has been decomposed by means of an Cholesky Decomposition.

The functions fail if ...

The first four functions report failure via the info argument, the fifth function throws a std::invalid_argument exception in case of an error.


Linear System Solver for Triangular Matrices

The following functions provide an interface for the LAPACK functions strsv(), dtrsv(), ctrsv(), and ztrsv():

namespace blaze {
void trsv( char uplo, char trans, char diag, int n, const float* A, int lda, float* x, int incX );
void trsv( char uplo, char trans, char diag, int n, const double* A, int lda, double* x, int incX );
void trsv( char uplo, char trans, char diag, int n, const complex<float>* A, int lda, complex<float>* x, int incX );
void trsv( char uplo, char trans, char diag, int n, const complex<double>* A, int lda, complex<double>* x, int incX );
template< typename MT, bool SO, typename VT, bool TF >
void trsv( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b, char uplo, char trans, char diag );
} // namespace blaze

If the function exits successfully, the vector b or the matrix B contain the solution(s) of the linear system of equations.

The functions fail if ...

The last function throws a std::invalid_argument exception in case of an error. Note that none of the functions does perform any test for singularity or near-singularity. Such tests must be performed prior to calling this function!


Previous: BLAS Functions     Next: Configuration Files