TRIQS/nda 1.3.0
Multi-dimensional array library for C++
Loading...
Searching...
No Matches
LAPACK interface

Detailed Description

Interface to parts of the LAPACK library.

Classes

class  nda::lapack::gelss_worker< T >
 Worker class for solving linear least square problems. More...
 
struct  nda::lapack::gelss_worker_hermitian
 Worker class for solving linear least square problems for hermitian tail-fitting. More...
 

Functions

template<MemoryMatrix A, MemoryArray B, MemoryVector S>
requires (have_same_value_type_v<A, B> and mem::on_host<A, B, S> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::gelss (A &&a, B &&b, S &&s, double rcond, int &rank)
 Interface to the LAPACK gelss routine.
 
template<MemoryMatrix A, MemoryVector JPVT, MemoryVector TAU>
requires (mem::on_host<A> and is_blas_lapack_v<get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, JPVT, TAU>)
int nda::lapack::geqp3 (A &&a, JPVT &&jpvt, TAU &&tau)
 Interface to the LAPACK geqp3 routine.
 
template<MemoryMatrix A, MemoryVector S, MemoryMatrix U, MemoryMatrix VT>
requires (have_same_value_type_v<A, U, VT> and mem::have_compatible_addr_space<A, S, U, VT> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::gesvd (A &&a, S &&s, U &&u, VT &&vt)
 Interface to the LAPACK gesvd routine.
 
template<MemoryMatrix A, MemoryVector IPIV>
requires (mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getrf (A &&a, IPIV &&ipiv)
 Interface to the LAPACK getrf routine.
 
template<MemoryMatrix A, MemoryVector IPIV>
requires (mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getri (A &&a, IPIV const &ipiv)
 Interface to the LAPACK getri routine.
 
template<MemoryMatrix A, MemoryMatrix B, MemoryVector IPIV>
requires (have_same_value_type_v<A, B> and mem::have_compatible_addr_space<A, B, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getrs (A const &a, B &&b, IPIV const &ipiv)
 Interface to the LAPACK getrs routine.
 
template<MemoryVector DL, MemoryVector D, MemoryVector DU, MemoryArray B>
requires (have_same_value_type_v<DL, D, DU, B> and mem::on_host<DL, D, DU, B> and is_blas_lapack_v<get_value_t<DL>>)
int nda::lapack::gtsv (DL &&dl, D &&d, DU &&du, B &&b)
 Interface to the LAPACK gtsv routine.
 
template<MemoryMatrix A, MemoryVector TAU>
requires (mem::on_host<A> and std::is_same_v<double, get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, TAU>)
int nda::lapack::orgqr (A &&a, TAU &&tau)
 Interface to the LAPACK orgqr routine.
 
template<MemoryMatrix A, MemoryVector TAU>
requires (mem::on_host<A> and std::is_same_v<std::complex<double>, get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, TAU>)
int nda::lapack::ungqr (A &&a, TAU &&tau)
 Interface to the LAPACK ungqr routine.
 

Function Documentation

◆ gelss()

template<MemoryMatrix A, MemoryArray B, MemoryVector S>
requires (have_same_value_type_v<A, B> and mem::on_host<A, B, S> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::gelss ( A && a,
B && b,
S && s,
double rcond,
int & rank )

#include <nda/lapack/gelss.hpp>

Interface to the LAPACK gelss routine.

Computes the minimum norm solution to a complex linear least squares problem:

\[ \min_x | \mathbf{b} - \mathbf{A x} |_2 \]

using the singular value decomposition (SVD) of \( \mathbf{A} \). \( \mathbf{A} \) is an m-by-n matrix which may be rank-deficient.

Several right hand side vectors \( \mathbf{b} \) and solution vectors \( \mathbf{x} \) can be handled in a single call; they are stored as the columns of the m-by-nrhs right hand side matrix \( \mathbf{B} \) and the n-by-nrhs solution matrix \( \mathbf{X} \).

The effective rank of \( \mathbf{A} \) is determined by treating as zero those singular values which are less than rcond times the largest singular value.

Template Parameters
Anda::MemoryMatrix type.
Bnda::MemoryArray type.
Snda::MemoryVector type.
Parameters
aInput/output matrix. On entry, the m-by-n matrix \( \mathbf{A} \). On exit, the first min(m,n) rows of \( \mathbf{A} \) are overwritten with its right singular vectors, stored rowwise.
bInput/output array. On entry, the m-by-nrhs right hand side matrix \( \mathbf{B} \). On exit, \( \mathbf{B} \) is overwritten by the n-by-nrhs solution matrix \( \mathbf{X} \). If m >= n and RANK == n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of the modulus of elements n+1:m in that column.
sOutput vector. The singular values of \( \mathbf{A} \) in decreasing order. The condition number of A in the 2-norm is s(1)/s(min(m,n)).
rcondIt is used to determine the effective rank of \( \mathbf{A} \). Singular values s(i) <= rcond * s(1) are treated as zero. If rcond < 0, machine precision is used instead.
rankOutput variable of the effective rank of \( \mathbf{A} \), i.e., the number of singular values which are greater than rcond * s(1).
Returns
Integer return code.

Definition at line 76 of file gelss.hpp.

◆ geqp3()

template<MemoryMatrix A, MemoryVector JPVT, MemoryVector TAU>
requires (mem::on_host<A> and is_blas_lapack_v<get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, JPVT, TAU>)
int nda::lapack::geqp3 ( A && a,
JPVT && jpvt,
TAU && tau )

#include <nda/lapack/geqp3.hpp>

Interface to the LAPACK geqp3 routine.

Computes a QR factorization with column pivoting of a matrix \( \mathbf{A} \):

\[ \mathbf{A P} = \mathbf{Q R} \]

using Level 3 BLAS.

Template Parameters
Anda::MemoryMatrix type.
JPVTnda::MemoryVector type.
TAUnda::MemoryVector type.
Parameters
aInput/output matrix. On entry, the m-by-n matrix \( \mathbf{A} \). On exit, the upper triangle of the array contains the min(m,n)-by-n upper trapezoidal matrix \( \mathbf{R} \); the elements below the diagonal, together with the array tau, represent the unitary matrix \( \mathbf{Q} \) as a product of min(m,n) elementary reflectors.
jpvtInput/output vector. On entry, if jpvt(j) != 0, the j-th column of \( \mathbf{A} \) is permuted to the front of \( \mathbf{A P} \) (a leading column); if jpvt(j) == 0, the j-th column of \( \mathbf{A} \) is a free column. On exit, if jpvt(j) == k, then the j-th column of \( \mathbf{A P} \) was the the k-th column of \( \mathbf{A} \).
tauOutput vector. The scalar factors of the elementary reflectors.
Returns
Integer return code from the LAPACK call.

Definition at line 67 of file geqp3.hpp.

◆ gesvd()

template<MemoryMatrix A, MemoryVector S, MemoryMatrix U, MemoryMatrix VT>
requires (have_same_value_type_v<A, U, VT> and mem::have_compatible_addr_space<A, S, U, VT> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::gesvd ( A && a,
S && s,
U && u,
VT && vt )

#include <nda/lapack/gesvd.hpp>

Interface to the LAPACK gesvd routine.

Computes the singular value decomposition (SVD) of a complex m-by-n matrix \( \mathbf{A} \), optionally computing the left and/or right singular vectors. The SVD is written

\[ \mathbf{A} = \mathbf{U} \mathbf{S} \mathbf{V}^H \]

where \( \mathbf{S} \) is an m-by-n matrix which is zero except for its min(m,n) diagonal elements, \( \mathbf{U} \) is an m-by-m unitary matrix, and \( \mathbf{V} \) is an n-by-n unitary matrix. The diagonal elements of \( \mathbf{S} \) are the singular values of \( \mathbf{A} \); they are real and non-negative, and are returned in descending order. The first min(m,n) columns of \( \mathbf{U} \) and \( \mathbf{V} \) are the left and right singular vectors of \( \mathbf{A} \).

Note that the routine returns \( \mathbf{V}^H \), not \( \mathbf{V} \).

Template Parameters
Anda::MemoryMatrix type.
Snda::MemoryVector type.
Unda::MemoryMatrix type.
VTnda::MemoryMatrix type.
Parameters
aInput/output matrix. On entry, the m-by-n matrix \( \mathbf{A} \). On exit, the contents of \( \mathbf{A} \) are destroyed.
sOutput vector. The singular values of \( \mathbf{A} \), sorted so that s(i) >= s(i+1).
uOutput matrix. It contains the m-by-m unitary matrix \( \mathbf{U} \).
vtOutput matrix. It contains contains the n-by-n unitary matrix \( \mathbf{V}^H \).
Returns
Integer return code from the LAPACK call.

Definition at line 75 of file gesvd.hpp.

◆ getrf()

template<MemoryMatrix A, MemoryVector IPIV>
requires (mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getrf ( A && a,
IPIV && ipiv )

#include <nda/lapack/getrf.hpp>

Interface to the LAPACK getrf routine.

Computes an LU factorization of a general m-by-n matrix \( \mathbf{A} \) using partial pivoting with row interchanges.

The factorization has the form

\[ \mathbf{A} = \mathbf{P L U} \]

where \( \mathbf{P} \) is a permutation matrix, \( \mathbf{L} \) is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and \( \mathbf{U} \) is upper triangular (upper trapezoidal if m < n).

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

Template Parameters
Anda::MemoryMatrix type.
IPIVnda::MemoryVector type.
Parameters
aInput/output matrix. On entry, the m-by-n matrix to be factored. On exit, the factors \( \mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \); the unit diagonal elements of \( \mathbf{L} \) are not stored.
ipivOutput vector. The pivot indices from getrf, i.e. for 1 <= i <= N, row i of the matrix was interchanged with row ipiv(i).
Returns
Integer return code from the LAPACK call.

Definition at line 62 of file getrf.hpp.

◆ getri()

template<MemoryMatrix A, MemoryVector IPIV>
requires (mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getri ( A && a,
IPIV const & ipiv )

#include <nda/lapack/getri.hpp>

Interface to the LAPACK getri routine.

Computes the inverse of a matrix using the LU factorization computed by getrf.

This method inverts \( \mathbf{U} \) and then computes \( \mathrm{inv}(\mathbf{A}) \) by solving the system \( \mathrm{inv}(\mathbf{A}) L = \mathrm{inv}(\mathbf{U}) \) for \( \mathrm{inv}(\mathbf{A}) \).

Template Parameters
Anda::MemoryMatrix type.
IPIVnda::MemoryVector type.
Parameters
aInput/output matrix. On entry, the factors \( \mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \) as computed by getrf. On exit, if INFO == 0, the inverse of the original matrix \( \mathbf{A} \).
ipivInput vector. The pivot indices from getrf, i.e. for 1 <= i <= N, row i of the matrix was interchanged with row ipiv(i).
Returns
Integer return code from the LAPACK call.

Definition at line 59 of file getri.hpp.

◆ getrs()

template<MemoryMatrix A, MemoryMatrix B, MemoryVector IPIV>
requires (have_same_value_type_v<A, B> and mem::have_compatible_addr_space<A, B, IPIV> and is_blas_lapack_v<get_value_t<A>>)
int nda::lapack::getrs ( A const & a,
B && b,
IPIV const & ipiv )

#include <nda/lapack/getrs.hpp>

Interface to the LAPACK getrs routine.

Solves a system of linear equations

  • \( \mathbf{A X} = \mathbf{B} \),
  • \( \mathbf{A}^T \mathbf{X} = \mathbf{B} \) or
  • \( \mathbf{A}^H \mathbf{X} = \mathbf{B} \)

with a general n-by-n matrix \( \mathbf{A} \) using the LU factorization computed by getrf.

Template Parameters
Anda::MemoryMatrix type.
Bnda::MemoryMatrix type.
Bnda::MemoryVector type.
Parameters
aInput matrix. The factors \( \mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \) as computed by getrf.
bInput/output matrix. On entry, the right hand side matrix \( \mathbf{B} \). On exit, the solution matrix \( \mathbf{X} \).
ipivInput vector. The pivot indices from getrf, i.e. for 1 <= i <= N, row i of the matrix was interchanged with row ipiv(i).
Returns
Integer return code from the LAPACK call.

Definition at line 64 of file getrs.hpp.

◆ gtsv()

template<MemoryVector DL, MemoryVector D, MemoryVector DU, MemoryArray B>
requires (have_same_value_type_v<DL, D, DU, B> and mem::on_host<DL, D, DU, B> and is_blas_lapack_v<get_value_t<DL>>)
int nda::lapack::gtsv ( DL && dl,
D && d,
DU && du,
B && b )

#include <nda/lapack/gtsv.hpp>

Interface to the LAPACK gtsv routine.

Solves the equation

\[ \mathbf{A} \mathbf{X} = \mathbf{B}, \]

where \( \mathbf{A} \) is an n-by-n tridiagonal matrix, by Gaussian elimination with partial pivoting.

Note that the equation \( \mathbf{A}^T \mathbf{X} = \mathbf{B} \) may be solved by interchanging the order of the arguments containing the subdiagonal elements.

Template Parameters
DLnda::MemoryVector type.
Dnda::MemoryVector type.
DUnda::MemoryVector type.
Bnda::MemoryArray type.
Parameters
dlInput/Output vector. On entry, it must contain the (n-1) subdiagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the (n-2) elements of the second superdiagonal of the upper triangular matrix \( \mathbf{U} \) from the LU factorization of \( \mathbf{A} \).
dInput/Output vector. On entry, it must contain the diagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the n diagonal elements of \( \mathbf{U} \).
duInput/Output vector. On entry, it must contain the (n-1) superdiagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the (n-1) elements of the first superdiagonal of \( \mathbf{U} \) .
bInput/Output array. On entry, the n-by-nrhs right hand side matrix \( \mathbf{B} \). On exit, if INFO == 0, the n-by-nrhs solution matrix \( \mathbf{X} \).
Returns
Integer return code from the LAPACK call.

Definition at line 62 of file gtsv.hpp.

◆ orgqr()

template<MemoryMatrix A, MemoryVector TAU>
requires (mem::on_host<A> and std::is_same_v<double, get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, TAU>)
int nda::lapack::orgqr ( A && a,
TAU && tau )

#include <nda/lapack/orgqr.hpp>

Interface to the LAPACK orgqr routine.

Generates an m-by-n real matrix \( \mathbf{Q} \) with orthonormal columns, which is defined as the first n columns of a product of k elementary reflectors of order m:

\[ \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; , \]

as returned by geqrf.

Template Parameters
Anda::MemoryMatrix with double value type.
TAUnda::MemoryVector with double value type.
Parameters
aInput/output matrix. On entry, the i-th column must contain the vector which defines the elementary reflector \( H(i) \; , i = 1,2,...,k \), as returned by geqrf in the first k columns. On exit, the m-by-n matrix \( \mathbf{Q} \).
tauInput vector. tau(i) must contain the scalar factor of the elementary reflector \( \mathbf{H}(i) \), as returned by geqrf.
Returns
Integer return code from the LAPACK call.

Definition at line 63 of file orgqr.hpp.

◆ ungqr()

template<MemoryMatrix A, MemoryVector TAU>
requires (mem::on_host<A> and std::is_same_v<std::complex<double>, get_value_t<A>> and have_same_value_type_v<A, TAU> and mem::have_compatible_addr_space<A, TAU>)
int nda::lapack::ungqr ( A && a,
TAU && tau )

#include <nda/lapack/ungqr.hpp>

Interface to the LAPACK ungqr routine.

Generates an m-by-n complex matrix \( \mathbf{Q} \) with orthonormal columns, which is defined as the first n columns of a product of k elementary reflectors of order m

\[ \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; , \]

as returned by geqrf.

Template Parameters
Anda::MemoryMatrix with complex value type.
TAUnda::MemoryVector with complex value type.
Parameters
aInput/output matrix. On entry, the i-th column must contain the vector which defines the elementary reflector \( H(i) \; , i = 1,2,...,K \), as returned by geqrf in the first k columns. On exit, the m-by-n matrix \( \mathbf{Q} \).
tauInput vector. tau(i) must contain the scalar factor of the elementary reflector \( \mathbf{H}(i) \), as returned by geqrf.
Returns
Integer return code from the LAPACK call.

Definition at line 63 of file ungqr.hpp.