<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
Wrappers to variants of the computational routines

Wrappers that receive a parameter to select a variant of the computational routine. More...

Functions

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::bidiag (matrix_t &A, vector_t &tauv, vector_t &tauw, const BidiagOpts &opts={})
 Reduces a general m by n matrix A to an upper real bidiagonal form B by a unitary transformation:
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::bidiag_work (matrix_t &A, vector_t &tauv, vector_t &tauw, work_t &work, const BidiagOpts &opts={})
 Reduces a general m by n matrix A to an upper real bidiagonal form B by a unitary transformation: Workspace is provided as an argument.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::gen_householder_q (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau, const GenHouseholderQOpts &opts={})
 Generates a matrix Q that is the product of elementary reflectors.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t>
int tlapack::gen_householder_q_work (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau, work_t &work, const GenHouseholderQOpts &opts={})
 Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR piv_t>
int tlapack::getrf (matrix_t &A, piv_t &piv, const GetrfOpts &opts={})
 getrf computes an LU factorization of a general m-by-n matrix A.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR piv_t>
int tlapack::getri (matrix_t &A, const piv_t &piv, const GetriOpts &opts={})
 getri computes inverse of a general n-by-n matrix A
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR piv_t, TLAPACK_WORKSPACE work_t>
int tlapack::getri_work (matrix_t &A, const piv_t &piv, work_t &work, const GetriOpts &opts={})
 getri computes inverse of a general n-by-n matrix A Workspace is provided as an argument.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::hessenberg (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau, const HessenbergOpts &opts={})
 Reduces a general square matrix to upper Hessenberg form.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::hessenberg_work (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau, work_t &work, const HessenbergOpts &opts={})
 Reduces a general square matrix to upper Hessenberg form. Workspace is provided as an argument.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_lq (matrix_t &A, vector_t &tau, const HouseholderLQOpts &opts={})
 Computes a LQ factorization of an m-by-n matrix A.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_lq_work (matrix_t &A, vector_t &tau, workspace_t &work, const HouseholderLQOpts &opts={})
 Computes a LQ factorization of an m-by-n matrix A. Workspace is provided as an argument.
 
template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::householder_q_mul (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const vector_t &tau, matrixC_t &C, const HouseholderQMulOpts &opts={})
 Applies unitary matrix Q to a matrix C.
 
template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t>
int tlapack::householder_q_mul_work (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const vector_t &tau, matrixC_t &C, work_t &work, const HouseholderQMulOpts &opts={})
 Applies unitary matrix Q to a matrix C. Workspace is provided as an argument.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_ql (matrix_t &A, vector_t &tau, const HouseholderQLOpts &opts={})
 Computes a QL factorization of an m-by-n matrix A.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::householder_ql_work (matrix_t &A, vector_t &tau, work_t &work, const HouseholderQLOpts &opts={})
 Computes a QL factorization of an m-by-n matrix A. Workspace is provided as an argument.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_qr (matrix_t &A, vector_t &tau, const HouseholderQROpts &opts={})
 Computes a QR factorization of an m-by-n matrix A.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_qr_work (matrix_t &A, vector_t &tau, workspace_t &work, const HouseholderQROpts &opts={})
 Computes a QR factorization of an m-by-n matrix A. Workspace is provided as an argument.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_rq (matrix_t &A, vector_t &tau, const HouseholderRQOpts &opts={})
 Computes a RQ factorization of an m-by-n matrix A.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_rq_work (matrix_t &A, vector_t &tau, workspace_t &work, const HouseholderRQOpts &opts={})
 Computes a RQ factorization of an m-by-n matrix A. Workspace is provided as an argument.
 
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t>
int tlapack::potrf (uplo_t uplo, matrix_t &A, const PotrfOpts &opts={})
 Computes the Cholesky factorization of a Hermitian positive definite matrix A.
 
template<class d_t , class e_t , enable_if_t< is_same_v< real_type< type_t< d_t > >, real_type< type_t< e_t > > >, int > = 0>
int tlapack::pttrf (d_t &D, e_t &E, const EcOpts &opts={})
 Computes the Cholesky factorization of a Hermitian positive definite tridiagonal matrix A.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::qr_iteration (bool want_t, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &w, matrix_t &Z, QRIterationOpts &opts)
 Computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix.
 
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::qr_iteration_work (bool want_t, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &w, matrix_t &Z, work_t &work, QRIterationOpts &opts)
 Computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix. Workspace is provided as an argument.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::unghr (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, const vector_t &tau)
 Generates a m-by-n matrix Q with orthogonal columns.
 

Detailed Description

Wrappers that receive a parameter to select a variant of the computational routine.


Function Documentation

◆ bidiag()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::bidiag ( matrix_t A,
vector_t tauv,
vector_t tauw,
const BidiagOpts opts = {} 
)

Reduces a general m by n matrix A to an upper real bidiagonal form B by a unitary transformation:

\[ Q**H * A * P = B. \]

The matrices Q and P are represented as products of elementary reflectors:

If m >= n,

\[ Q = H(1) H(2) . . . H(n) and P = G(1) G(2) . . . G(n-1) \]

Each H(i) and G(i) has the form:

\[ H(j) = I - tauv * v * v**H and G(j) = I - tauw * w * w**H \]

where tauv and tauw are scalars, and v and w are vectors; v(1:j-1) = 0, v(j) = 1, and v(j+1:m) is stored on exit in A(j+1:m,j); w(1:j) = 0, w(j+1) = 1, and w(j+2:n) is stored on exit in A(j,i+2:n); tauv is stored in tauv(j) and tauw in tauw(j).

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On entry, the m by n general matrix to be reduced. On exit,
  • if m >= n, the diagonal and the first superdiagonal are overwritten with the upper bidiagonal matrix B; the elements below the diagonal, with the array tauv, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the first superdiagonal, with the array tauw, represent the unitary matrix P as a product of elementary reflectors.
  • if m < n, the diagonal and the first superdiagonal are overwritten with the lower bidiagonal matrix B; the elements below the first subdiagonal, with the array tauv, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the diagonal, with the array tauw, represent the unitary matrix P as a product of elementary reflectors.
[out]tauvvector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Q.
[out]tauwvector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix P.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.

◆ bidiag_work()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::bidiag_work ( matrix_t A,
vector_t tauv,
vector_t tauw,
work_t work,
const BidiagOpts opts = {} 
)

Reduces a general m by n matrix A to an upper real bidiagonal form B by a unitary transformation: Workspace is provided as an argument.

\[ Q**H * A * P = B. \]

The matrices Q and P are represented as products of elementary reflectors:

If m >= n,

\[ Q = H(1) H(2) . . . H(n) and P = G(1) G(2) . . . G(n-1) \]

Each H(i) and G(i) has the form:

\[ H(j) = I - tauv * v * v**H and G(j) = I - tauw * w * w**H \]

where tauv and tauw are scalars, and v and w are vectors; v(1:j-1) = 0, v(j) = 1, and v(j+1:m) is stored on exit in A(j+1:m,j); w(1:j) = 0, w(j+1) = 1, and w(j+2:n) is stored on exit in A(j,i+2:n); tauv is stored in tauv(j) and tauw in tauw(j).

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On entry, the m by n general matrix to be reduced. On exit,
  • if m >= n, the diagonal and the first superdiagonal are overwritten with the upper bidiagonal matrix B; the elements below the diagonal, with the array tauv, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the first superdiagonal, with the array tauw, represent the unitary matrix P as a product of elementary reflectors.
  • if m < n, the diagonal and the first superdiagonal are overwritten with the lower bidiagonal matrix B; the elements below the first subdiagonal, with the array tauv, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the diagonal, with the array tauw, represent the unitary matrix P as a product of elementary reflectors.
[out]tauvvector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Q.
[out]tauwvector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix P.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
workWorkspace. Use the workspace query to determine the size needed.

◆ gen_householder_q()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::gen_householder_q ( direction_t  direction,
storage_t  storeMode,
matrix_t A,
const vector_t tau,
const GenHouseholderQOpts opts = {} 
)

Generates a matrix Q that is the product of elementary reflectors.

Parameters
[in]directionIndicates how Q is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in,out]Am-by-n matrix. On entry,
[in]tauVector of length k. Scalar factors of the elementary reflectors.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
Returns
0 if success.
Further Details

The shape of the matrix V and the storage of the vectors which define the \(H_i\) is best illustrated by the following example with k = 3. The elements equal to 1 are not accessed. The rest of the matrix is not used.

direction = Forward and          direction = Forward and
storeMode = Columnwise:             storeMode = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storeMode = Columnwise:             storeMode = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )

◆ gen_householder_q_work()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t>
int tlapack::gen_householder_q_work ( direction_t  direction,
storage_t  storeMode,
matrix_t A,
const vector_t tau,
work_t work,
const GenHouseholderQOpts opts = {} 
)

Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument.

Parameters
[in]directionIndicates how Q is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in,out]Am-by-n matrix. On entry,
[in]tauVector of length k. Scalar factors of the elementary reflectors.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
Returns
0 if success.
Further Details

The shape of the matrix V and the storage of the vectors which define the \(H_i\) is best illustrated by the following example with k = 3. The elements equal to 1 are not accessed. The rest of the matrix is not used.

direction = Forward and          direction = Forward and
storeMode = Columnwise:             storeMode = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storeMode = Columnwise:             storeMode = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )
Parameters
workWorkspace. Use the workspace query to determine the size needed.

◆ getrf()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR piv_t>
int tlapack::getrf ( matrix_t A,
piv_t piv,
const GetrfOpts opts = {} 
)

getrf computes an LU factorization of a general m-by-n matrix A.

The factorization has the form

\[ P A = L U \]

where P is a permutation matrix constructed from our piv vector, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

Returns
0 if success
i+1 if failed to compute the LU on iteration i
Parameters
[in,out]Am-by-n matrix. On exit, the factors L and U from the factorization A=PLU; the unit diagonal elements of L are not stored.
[in,out]pivis a k-by-1 integer vector where k=min(m,n) and piv[i]=j where i<=j<=k-1, which means in the i-th iteration of the algorithm, the j-th row needs to be swapped with i
[in]optsOptions.
  • variant:
    • Recursive = 'R',
    • Level0 = '0'
Note
To construct L and U, one proceeds as in the following steps
  1. Set matrices L m-by-k, and U k-by-n be to matrices with all zeros, where k=min(m,n)
  2. Set elements on the diagonal of L to 1
  3. below the diagonal of A will be copied to L
  4. On and above the diagonal of A will be copied to U

◆ getri()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR piv_t>
int tlapack::getri ( matrix_t A,
const piv_t piv,
const GetriOpts opts = {} 
)

getri computes inverse of a general n-by-n matrix A

Returns
= 0: successful exit
= i+1: if U(i,i) is exactly zero. The triangular matrix is singular and its inverse can not be computed.
Parameters
[in,out]An-by-n matrix. On entry, the factors L and U from the factorization P A = L U. L is stored in the lower triangle of A, the unit diagonal elements of L are not stored. U is stored in the upper triangle of A. On exit, inverse of A is overwritten on A.
[in]pivpivot vector of size at least n.
[in]optsOptions.
  • opts.variant:
    • UILI = 'D', ///< Method D from doi:10.1137/1.9780898718027
    • UXLI = 'C' ///< Method C from doi:10.1137/1.9780898718027

◆ getri_work()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR piv_t, TLAPACK_WORKSPACE work_t>
int tlapack::getri_work ( matrix_t A,
const piv_t piv,
work_t work,
const GetriOpts opts = {} 
)

getri computes inverse of a general n-by-n matrix A Workspace is provided as an argument.

Returns
= 0: successful exit
= i+1: if U(i,i) is exactly zero. The triangular matrix is singular and its inverse can not be computed.
Parameters
[in,out]An-by-n matrix. On entry, the factors L and U from the factorization P A = L U. L is stored in the lower triangle of A, the unit diagonal elements of L are not stored. U is stored in the upper triangle of A. On exit, inverse of A is overwritten on A.
[in]pivpivot vector of size at least n.
[in]optsOptions.
  • opts.variant:
    • UILI = 'D', ///< Method D from doi:10.1137/1.9780898718027
    • UXLI = 'C' ///< Method C from doi:10.1137/1.9780898718027
workWorkspace. Use the workspace query to determine the size needed.

◆ hessenberg()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::hessenberg ( size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t tau,
const HessenbergOpts opts = {} 
)

Reduces a general square matrix to upper Hessenberg form.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_ilo H_ilo+1 ... H_ihi, \]

Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i] = 0; v[i+1] = 1, \]

with v[i+2] through v[ihi] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Returns
0 if success
Parameters
[in]ilointeger
[in]ihiinteger It is assumed that A is already upper Hessenberg in columns 0:ilo and rows ihi:n and is already upper triangular in columns ihi+1:n and rows 0:ilo. 0 <= ilo <= ihi <= max(1,n).
[in,out]An-by-n matrix. On entry, the n by n general matrix to be reduced. On exit, the upper triangle and the first subdiagonal of A are overwritten with the upper Hessenberg matrix H, and the elements below the first subdiagonal, with the array TAU, represent the orthogonal matrix Q as a product of elementary reflectors. See Further Details.
[out]tauReal vector of length n-1. The scalar factors of the elementary reflectors.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.

◆ hessenberg_work()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::hessenberg_work ( size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t tau,
work_t work,
const HessenbergOpts opts = {} 
)

Reduces a general square matrix to upper Hessenberg form. Workspace is provided as an argument.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_ilo H_ilo+1 ... H_ihi, \]

Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i] = 0; v[i+1] = 1, \]

with v[i+2] through v[ihi] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Returns
0 if success
Parameters
[in]ilointeger
[in]ihiinteger It is assumed that A is already upper Hessenberg in columns 0:ilo and rows ihi:n and is already upper triangular in columns ihi+1:n and rows 0:ilo. 0 <= ilo <= ihi <= max(1,n).
[in,out]An-by-n matrix. On entry, the n by n general matrix to be reduced. On exit, the upper triangle and the first subdiagonal of A are overwritten with the upper Hessenberg matrix H, and the elements below the first subdiagonal, with the array TAU, represent the orthogonal matrix Q as a product of elementary reflectors. See Further Details.
[out]tauReal vector of length n-1. The scalar factors of the elementary reflectors.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
workWorkspace. Use the workspace query to determine the size needed.

◆ householder_lq()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_lq ( matrix_t A,
vector_t tau,
const HouseholderLQOpts opts = {} 
)

Computes a LQ factorization of an m-by-n matrix A.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with conj(v[i+1]) through conj(v[n-1]) stored on exit above the diagonal in the ith row of A, and conj(tau) in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, the elements on and below the diagonal of the array contain the m-by-k lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ householder_lq_work()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_lq_work ( matrix_t A,
vector_t tau,
workspace_t work,
const HouseholderLQOpts opts = {} 
)

Computes a LQ factorization of an m-by-n matrix A. Workspace is provided as an argument.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with conj(v[i+1]) through conj(v[n-1]) stored on exit above the diagonal in the ith row of A, and conj(tau) in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, the elements on and below the diagonal of the array contain the m-by-k lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.
workWorkspace. Use the workspace query to determine the size needed.

◆ householder_q_mul()

template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::householder_q_mul ( side_t  side,
trans_t  trans,
direction_t  direction,
storage_t  storeMode,
const matrixV_t V,
const vector_t tau,
matrixC_t C,
const HouseholderQMulOpts opts = {} 
)

Applies unitary matrix Q to a matrix C.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with v[i+1] through v[m-1] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]directionIndicates how Q is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]V
[in]tauVector of length k. Scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
Returns
0 if success.
Further Details

The shape of the matrix V and the storage of the vectors which define the \(H_i\) is best illustrated by the following example with k = 3. The elements equal to 1 are not accessed. The rest of the matrix is not used.

direction = Forward and          direction = Forward and
storeMode = Columnwise:             storeMode = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storeMode = Columnwise:             storeMode = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )

◆ householder_q_mul_work()

template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t>
int tlapack::householder_q_mul_work ( side_t  side,
trans_t  trans,
direction_t  direction,
storage_t  storeMode,
const matrixV_t V,
const vector_t tau,
matrixC_t C,
work_t work,
const HouseholderQMulOpts opts = {} 
)

Applies unitary matrix Q to a matrix C. Workspace is provided as an argument.

- side = Side::Left & trans = Op::NoTrans: \(C := Q C\);

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with v[i+1] through v[m-1] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]directionIndicates how Q is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]V
[in]tauVector of length k. Scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
Returns
0 if success.
Further Details

The shape of the matrix V and the storage of the vectors which define the \(H_i\) is best illustrated by the following example with k = 3. The elements equal to 1 are not accessed. The rest of the matrix is not used.

direction = Forward and          direction = Forward and
storeMode = Columnwise:             storeMode = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storeMode = Columnwise:             storeMode = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )
Parameters
workWorkspace. Use the workspace query to determine the size needed.

◆ householder_ql()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_ql ( matrix_t A,
vector_t tau,
const HouseholderQLOpts opts = {} 
)

Computes a QL factorization of an m-by-n matrix A.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_k ... H_2 H_1, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[m-k+i+1:m] = 0; v[m-k+i-1] = 1, \]

with v[1] through v[m-k+i-1] stored on exit in A(0:m-k+i-1,n-k+i), and tau in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, if m >= n, the lower triangle of A(m-n:m,0:n) contains the n by n lower triangular matrix L; If m <= n, the elements on and below the (n-m)-th superdiagonal contain the m by n lower trapezoidal matrix L the remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ householder_ql_work()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t>
int tlapack::householder_ql_work ( matrix_t A,
vector_t tau,
work_t work,
const HouseholderQLOpts opts = {} 
)

Computes a QL factorization of an m-by-n matrix A. Workspace is provided as an argument.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_k ... H_2 H_1, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[m-k+i+1:m] = 0; v[m-k+i-1] = 1, \]

with v[1] through v[m-k+i-1] stored on exit in A(0:m-k+i-1,n-k+i), and tau in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, if m >= n, the lower triangle of A(m-n:m,0:n) contains the n by n lower triangular matrix L; If m <= n, the elements on and below the (n-m)-th superdiagonal contain the m by n lower trapezoidal matrix L the remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.
workWorkspace. Use the workspace query to determine the size needed.

◆ householder_qr()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_qr ( matrix_t A,
vector_t tau,
const HouseholderQROpts opts = {} 
)

Computes a QR factorization of an m-by-n matrix A.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with v[i+1] through v[m-1] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, the elements on and above the diagonal of the array contain the k-by-n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ householder_qr_work()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_qr_work ( matrix_t A,
vector_t tau,
workspace_t work,
const HouseholderQROpts opts = {} 
)

Computes a QR factorization of an m-by-n matrix A. Workspace is provided as an argument.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_1 H_2 ... H_k, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[0] = v[1] = ... = v[i-1] = 0; v[i] = 1, \]

with v[i+1] through v[m-1] stored on exit below the diagonal in the ith column of A, and tau in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, the elements on and above the diagonal of the array contain the k-by-n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.
workWorkspace. Use the workspace query to determine the size needed.

◆ householder_rq()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::householder_rq ( matrix_t A,
vector_t tau,
const HouseholderRQOpts opts = {} 
)

Computes a RQ factorization of an m-by-n matrix A.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_k ... H_2 H_1, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[n-k+i+1:n] = 0; v[n-k+i-1] = 1, \]

with conj(v[1]) through conj(v[n-k+i-1]) stored on exit in A(m-k+i,0:n-k+i-1), and conj(tau) in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, if m <= n, the upper triangle of the subarray A(0:m,n-m:n) contains the m by m upper triangular matrix R; if m >= n, the elements on and above the (m-n)-th subdiagonal contain the m by n upper trapezoidal matrix R; the remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ householder_rq_work()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE workspace_t>
int tlapack::householder_rq_work ( matrix_t A,
vector_t tau,
workspace_t work,
const HouseholderRQOpts opts = {} 
)

Computes a RQ factorization of an m-by-n matrix A. Workspace is provided as an argument.

The matrix Q is represented as a product of elementary reflectors

\[ Q = H_k ... H_2 H_1, \]

where k = min(m,n). Each H_i has the form

\[ H_i = I - tau * v * v', \]

where tau is a scalar, and v is a vector with

\[ v[n-k+i+1:n] = 0; v[n-k+i-1] = 1, \]

with conj(v[1]) through conj(v[n-k+i-1]) stored on exit in A(m-k+i,0:n-k+i-1), and conj(tau) in tau[i].

Returns
0 if success
Parameters
[in,out]Am-by-n matrix. On exit, if m <= n, the upper triangle of the subarray A(0:m,n-m:n) contains the m by m upper triangular matrix R; if m >= n, the elements on and above the (m-n)-th subdiagonal contain the m by n upper trapezoidal matrix R; the remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length k. The scalar factors of the elementary reflectors.
[in]optsOptions.
workWorkspace. Use the workspace query to determine the size needed.

◆ potrf()

template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t>
int tlapack::potrf ( uplo_t  uplo,
matrix_t A,
const PotrfOpts opts = {} 
)

Computes the Cholesky factorization of a Hermitian positive definite matrix A.

The factorization has the form \(A = U^H U,\) if uplo = Upper, or \(A = L L^H,\) if uplo = Lower, where U is an upper triangular matrix and L is lower triangular.

Template Parameters
uplo_tAccess type: Upper or Lower. Either Uplo or any class that implements operator Uplo().
Parameters
[in]uplo
[in,out]AOn entry, the Hermitian matrix A of size n-by-n.
  • If uplo = Uplo::Upper, the strictly lower triangular part of A is not referenced.
  • If uplo = Uplo::Lower, the strictly upper triangular part of A is not referenced.
  • On successful exit, the factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H.\)
Parameters
[in]optsOptions. Define the behavior of checks for NaNs, and nb for potrf_blocked.
  • variant:
    • Recursive = 'R',
    • Blocked = 'B'
Returns
0: successful exit.
i, 0 < i <= n, if the leading minor of order i is not positive definite, and the factorization could not be completed.

◆ pttrf()

template<class d_t , class e_t , enable_if_t< is_same_v< real_type< type_t< d_t > >, real_type< type_t< e_t > > >, int > = 0>
int tlapack::pttrf ( d_t D,
e_t E,
const EcOpts opts = {} 
)

Computes the Cholesky factorization of a Hermitian positive definite tridiagonal matrix A.

The factorization has the form \(A = L D L^H\), where L is unit lower bidiagonal and D is diagonal.

Parameters
[in,out]DOn entry, the diagonal of \(A\).
  • On successful exit, the factor \(D\).
Parameters
[in,out]EOn entry, the first subdiagonal of \(A\).
  • On successful exit, the first subdiagonal of the factor \(L\).
Parameters
[in]optsOptions. Define the behavior of checks for NaNs.
Returns
0: successful exit.
i, 0 < i < n, if the leading minor of order i is not positive definite, and the factorization could not be completed.
n, if the leading minor of order n is not positive definite.

◆ qr_iteration()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::qr_iteration ( bool  want_t,
bool  want_z,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t w,
matrix_t Z,
QRIterationOpts opts 
)

Computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix.

Returns
0 if success
i if the QR algorithm failed to compute all the eigenvalues in a total of 30 iterations per eigenvalue. elements i:ihi of w contain those eigenvalues which have been successfully computed.
Parameters
[in]want_tbool. If true, the full Schur factor T will be computed.
[in]want_zbool. If true, the Schur vectors Z will be computed.
[in]ilointeger. Either ilo=0 or A(ilo,ilo-1) = 0.
[in]ihiinteger. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n.
[in,out]An by n matrix. On entry, the matrix A. On exit, if info=0 and want_t=true, the Schur factor T. T is quasi-triangular in rows and columns ilo:ihi, with the diagonal (block) entries in standard form (see above).
[out]wsize n vector. On exit, if info=0, w(ilo:ihi) contains the eigenvalues of A(ilo:ihi,ilo:ihi). The eigenvalues appear in the same order as the diagonal (block) entries of T.
[in,out]Zn by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A are accumulated into Z.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.

◆ qr_iteration_work()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::qr_iteration_work ( bool  want_t,
bool  want_z,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t w,
matrix_t Z,
work_t work,
QRIterationOpts opts 
)

Computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix. Workspace is provided as an argument.

Returns
0 if success
i if the QR algorithm failed to compute all the eigenvalues in a total of 30 iterations per eigenvalue. elements i:ihi of w contain those eigenvalues which have been successfully computed.
Parameters
[in]want_tbool. If true, the full Schur factor T will be computed.
[in]want_zbool. If true, the Schur vectors Z will be computed.
[in]ilointeger. Either ilo=0 or A(ilo,ilo-1) = 0.
[in]ihiinteger. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n.
[in,out]An by n matrix. On entry, the matrix A. On exit, if info=0 and want_t=true, the Schur factor T. T is quasi-triangular in rows and columns ilo:ihi, with the diagonal (block) entries in standard form (see above).
[out]wsize n vector. On exit, if info=0, w(ilo:ihi) contains the eigenvalues of A(ilo:ihi,ilo:ihi). The eigenvalues appear in the same order as the diagonal (block) entries of T.
[in,out]Zn by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A are accumulated into Z.
[in]optsOptions.
  • opts.variant: Variant of the algorithm to use.
workWorkspace. Use the workspace query to determine the size needed.

◆ unghr()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::unghr ( size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
const vector_t tau 
)

Generates a m-by-n matrix Q with orthogonal columns.

Parameters
[in]ilointeger
[in]ihiinteger ilo and ihi must have the same values as in the previous call to gehrd. Q is equal to the unit matrix except in the submatrix Q(ilo+1:ihi,ilo+1:ihi). 0 <= ilo <= ihi <= max(1,n).
[in,out]Am-by-n matrix. On entry, the vectors which define the elementary reflectors. On exit, the m-by-n matrix Q.
[in]tauReal vector of length n-1. The scalar factors of the elementary reflectors.