<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
Computational routines with automatic workspace allocation

Wrappers that allocate workspace and call the corresponding computational routine. More...

Functions

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
void tlapack::aggressive_early_deflation (bool want_t, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, size_type< matrix_t > nw, matrix_t &A, vector_t &s, matrix_t &Z, size_type< matrix_t > &ns, size_type< matrix_t > &nd, FrancisOpts &opts)
 aggressive_early_deflation accepts as input an upper Hessenberg matrix H and performs an orthogonal similarity transformation designed to detect and deflate fully converged eigenvalues from a trailing principal submatrix.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gebd2 (matrix_t &A, vector_t &tauv, vector_t &tauw)
 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>
int tlapack::gebrd (matrix_t &A, vector_t &tauv, vector_t &tauw, const GebrdOpts &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_VECTOR vector_t>
int tlapack::gehd2 (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau)
 Reduces a general square matrix to upper Hessenberg form.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::gehrd (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau, const GehrdOpts &opts={})
 Reduces a general square matrix to upper Hessenberg form.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gelq2 (matrix_t &A, vector_t &tauw)
 Computes an LQ factorization of a complex m-by-n matrix A using an unblocked algorithm.
 
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::gelqf (A_t &A, tau_t &tau, const GelqfOpts &opts={})
 Computes an LQ factorization of an m-by-n matrix A using a blocked algorithm.
 
template<TLAPACK_SMATRIX matrix_t>
int tlapack::gelqt (matrix_t &A, matrix_t &TT)
 Computes an LQ factorization of a complex m-by-n matrix A using a blocked algorithm.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::geql2 (matrix_t &A, vector_t &tau)
 Computes a QL factorization of a matrix A.
 
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::geqlf (A_t &A, tau_t &tau, const GeqlfOpts &opts={})
 Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::geqr2 (matrix_t &A, vector_t &tau)
 Computes a QR factorization of a matrix A.
 
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::geqrf (A_t &A, tau_t &tau, const GeqrfOpts &opts={})
 Computes a QR factorization of an m-by-n matrix A using a blocked algorithm.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gerq2 (matrix_t &A, vector_t &tau)
 Computes an RQ factorization of a matrix A.
 
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::gerqf (A_t &A, tau_t &tau, const GerqfOpts &opts={})
 Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm.
 
template<TLAPACK_SMATRIX matrix_t>
int tlapack::getri_uxli (matrix_t &A)
 getri computes inverse of a general n-by-n matrix A by solving for X in the following equation
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::multishift_qr (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, FrancisOpts &opts)
 multishift_qr computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix, using the multishift implicit QR algorithm with AED.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, bool > = true>
void tlapack::multishift_QR_sweep (bool want_t, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, const vector_t &s, matrix_t &Z)
 multishift_QR_sweep performs a single small-bulge multi-shift QR sweep.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR alpha_t, TLAPACK_VECTOR beta_t, enable_if_t< is_complex< type_t< alpha_t > >, bool > = true>
void tlapack::multishift_QZ_sweep (bool want_t, bool want_q, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, matrix_t &B, const alpha_t &alpha, const beta_t &beta, matrix_t &Q, matrix_t &Z)
 multishift_QR_sweep performs a single small-bulge multi-shift QR sweep.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ung2l (matrix_t &A, const vector_t &tau)
 Generates an m-by-n matrix Q with orthonormal columns, which is defined as the last n columns of a product of k elementary reflectors of order m.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ung2r (matrix_t &A, const vector_t &tau)
 Generates a matrix Q with orthogonal columns.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::ungbr_p (const size_type< matrix_t > k, matrix_t &A, const vector_t &tau, const UngbrOpts &opts={})
 Generates the unitary matrix P**H determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::ungbr_q (const size_type< matrix_t > k, matrix_t &A, const vector_t &tau, const UngbrOpts &opts={})
 Generates the unitary matrix Q determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ungl2 (matrix_t &Q, const vector_t &tauw)
 Generates all or part of the unitary matrix Q from an LQ factorization determined by gelq2 (unblocked algorithm).
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::ungq (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau, const UngqOpts &opts={})
 Generates a matrix Q that is the product of elementary reflectors.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::ungq_level2 (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau)
 Generates a matrix Q that is the product of elementary reflectors.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ungr2 (matrix_t &A, const vector_t &tau)
 Generates an m-by-n matrix Q with orthonormal columns, which is defined as the last m rows of a product of k elementary reflectors of order n.
 
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unm2l (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C)
 Applies unitary matrix Q from an QL factorization to a matrix C.
 
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unm2r (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C)
 Applies unitary matrix Q to a matrix C.
 
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::unmhr (Side side, Op trans, size_type< matrix_t > ilo, size_type< matrix_t > ihi, const matrix_t &A, const vector_t &tau, matrix_t &C)
 Applies unitary matrix Q to a matrix C.
 
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unml2 (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C)
 Applies unitary matrix Q from an LQ factorization 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>
int tlapack::unmq (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const vector_t &tau, matrixC_t &C, const UnmqOpts &opts={})
 Applies unitary matrix Q to a matrix C.
 
template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::unmq_level2 (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const vector_t &tau, matrixC_t &C)
 Applies unitary matrix Q to a matrix C.
 
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unmr2 (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C)
 Applies unitary matrix Q from an RQ factorization to a matrix C.
 

Detailed Description

Wrappers that allocate workspace and call the corresponding computational routine.


Function Documentation

◆ aggressive_early_deflation()

template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
void tlapack::aggressive_early_deflation ( bool  want_t,
bool  want_z,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
size_type< matrix_t nw,
matrix_t A,
vector_t s,
matrix_t Z,
size_type< matrix_t > &  ns,
size_type< matrix_t > &  nd,
FrancisOpts opts 
)

aggressive_early_deflation accepts as input an upper Hessenberg matrix H and performs an orthogonal similarity transformation designed to detect and deflate fully converged eigenvalues from a trailing principal submatrix.

On output H has been over- written by a new Hessenberg matrix that is a perturbation of an orthogonal similarity transformation of H. It is to be hoped that the final version of H has many zero subdiagonal entries.

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. ilo and ihi determine an isolated block in A.
[in]nwinteger. Desired window size to perform aggressive early deflation on. If the matrix is not large enough to provide the scratch space or if the isolated block is small, a smaller value may be used.
[in,out]An by n matrix. Hessenberg matrix on which AED will be performed
[out]ssize n vector. On exit, the entries s[ihi-nd-ns:ihi-nd] contain the unconverged eigenvalues that can be used a shifts. The entries s[ihi-nd:ihi] contain the converged eigenvalues. Entries outside the range s[ihi-nw:ihi] are not changed. The converged shifts are stored in the same positions as their correspinding diagonal elements in A.
[in,out]Zn by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A accumulated into Z.
[out]nsinteger. Number of eigenvalues available as shifts in s.
[out]ndinteger. Number of converged eigenvalues available as shifts in s.
[in,out]optsOptions.
  • Output parameters opts.n_aed, opts.n_sweep and opts.n_shifts_total are updated by the internal call to multishift_qr.

◆ gebd2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gebd2 ( matrix_t A,
vector_t tauv,
vector_t tauw 
)

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

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

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

If m >= n,

\[ Q = H(1) H(2) . . . H(n) and Z = 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 complex scalars, and v and w are complex 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.
  • 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 Z 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 Z as a product of elementary reflectors.
[out]tauvReal vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Q.
[out]tauwReal vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Z.

◆ gebrd()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::gebrd ( matrix_t A,
vector_t tauv,
vector_t tauw,
const GebrdOpts 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, \]

where m >= n.

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.

◆ gehd2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gehd2 ( size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t tau 
)

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.

◆ gehrd()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::gehrd ( size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
vector_t tau,
const GehrdOpts 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.

◆ gelq2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gelq2 ( matrix_t A,
vector_t tauw 
)

Computes an LQ factorization of a complex m-by-n matrix A using an unblocked algorithm.

The matrix Q is represented as a product of elementary reflectors.

\[ Q = H(k)**H ... H(2)**H H(1)**H, \]

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

\[ H(j) = I - tauw * w * w**H \]

where tauw is a complex scalar, and w is a complex vector with

\[ w[0] = w[1] = ... = w[j-1] = 0; w[j] = 1, \]

with w[j+1]**H through w[n]**H is stored on exit in the jth row of A, and tauw in tauw[j].

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 min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tauw, represent the unitary matrix Q as a product of elementary reflectors.
[out]tauwComplex vector of length min(m,n). The scalar factors of the elementary reflectors.

◆ gelqf()

template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::gelqf ( A_t A,
tau_t tau,
const GelqfOpts opts = {} 
)

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

The matrix Q is represented as a product of elementary reflectors.

\[ Q = H(k)**H ... H(2)**H H(1)**H, \]

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

\[ H(j) = I - tauw * w * w**H \]

where tauw is a scalar, and w is a vector with

\[ w[0] = w[1] = ... = w[j-1] = 0; w[j] = 1, \]

where w[j+1]**H through w[n]**H are stored on exit in the jth row of A.

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 min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tauw, represent the unitary matrix Q as a product of elementary reflectors.
[out]taumin(n,m) vector. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ gelqt()

template<TLAPACK_SMATRIX matrix_t>
int tlapack::gelqt ( matrix_t A,
matrix_t TT 
)

Computes an LQ factorization of a complex m-by-n matrix A using a blocked algorithm.

Stores the triangular factors for later use.

The matrix Q is represented as a product of elementary reflectors.

\[ Q = H(k)**H ... H(2)**H H(1)**H, \]

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

\[ H(j) = I - tauw * w * w**H \]

where tauw is a complex scalar, and w is a complex vector with

\[ w[0] = w[1] = ... = w[j-1] = 0; w[j] = 1, \]

with w[j+1]**H through w[n]**H is stored on exit in the jth row of A. tauw is stored in TT(j,i), where 0 <= i < nb and i = j (mod nb).

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 min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tauw, represent the unitary matrix Q as a product of elementary reflectors.
[out]TTm-by-nb matrix. In the representation of the block reflector. tauw[j] is stored in TT(j,i), where 0 <= i < nb and i = j (mod nb). On exit, TT( 0:k, 0:nb ) contains blocks used to build Q :

\[ Q^H = [ I - W(0:nb,0:k)^T * TT(0:nb,0:nb) * conj(W(0:nb,0:k)) ] * [ I - W(nb:2nb,0:k)^T * TT(nb:2nb,0:nb) * conj(W(nb:2nb,0:k)) ] * ... \]

For a good default of nb, see GelqfOpts

◆ geql2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::geql2 ( matrix_t A,
vector_t tau 
)

Computes a QL factorization of a 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[n-k+i-1] stored on exit below the diagonal 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 entry, the m by n matrix A. 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 min(m,n). The scalar factors of the elementary reflectors.

◆ geqlf()

template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::geqlf ( A_t A,
tau_t tau,
const GeqlfOpts opts = {} 
)

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

The matrix Q is represented as a product of elementary reflectors.

\[ Q = H(k)**H ... H(2)**H H(1)**H, \]

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

\[ H(j) = I - tauw * w * w**H \]

where tauw is a scalar, and w is a vector with

\[ w[0] = w[1] = ... = w[j-1] = 0; w[j] = 1, \]

where w[j+1]**H through w[n]**H are stored on exit in the jth row of A.

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 min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array tauw, represent the unitary matrix Q as a product of elementary reflectors.
[out]taumin(n,m) vector. The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ geqr2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::geqr2 ( matrix_t A,
vector_t tau 
)

Computes a QR factorization of a 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 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 unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length min(m,n). The scalar factors of the elementary reflectors.

◆ geqrf()

template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::geqrf ( A_t A,
tau_t tau,
const GeqrfOpts opts = {} 
)

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

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 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 unitary matrix Q as a product of elementary reflectors.
[out]tauReal vector of length min(m,n). The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ gerq2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::gerq2 ( matrix_t A,
vector_t tau 
)

Computes an RQ factorization of a 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[n-k+i+1:n] = 0; v[n-k+i-1] = 1, \]

with v[1] through v[n-k+i-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 entry, the m by n matrix A. 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 min(m,n). The scalar factors of the elementary reflectors.

◆ gerqf()

template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t>
int tlapack::gerqf ( A_t A,
tau_t tau,
const GerqfOpts opts = {} 
)

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

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[n-k+i+1:n] = 0; v[n-k+i-1] = 1, \]

with v[1] through v[n-k+i-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 entry, the m by n matrix A. 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 min(m,n). The scalar factors of the elementary reflectors.
[in]optsOptions.

◆ getri_uxli()

template<TLAPACK_SMATRIX matrix_t>
int tlapack::getri_uxli ( matrix_t A)

getri computes inverse of a general n-by-n matrix A by solving for X in the following equation

\[ U X L = I \]

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; unit diagonal is not stored. U is stored in the upper triangle of A. On exit, inverse of A is overwritten on A.

◆ multishift_qr()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0>
int tlapack::multishift_qr ( 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,
FrancisOpts opts 
)

multishift_qr computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix, using the multishift implicit QR algorithm with AED.

The Schur factorization is returned in standard form. For complex matrices this means that the matrix T is upper-triangular. The diagonal entries of T are also its eigenvalues. For real matrices, this means that the matrix T is block-triangular, with real eigenvalues appearing as 1x1 blocks on the diagonal and imaginary eigenvalues appearing as 2x2 blocks on the diagonal. All 2x2 blocks are normalized so that the diagonal entries are equal to the real part of the eigenvalue.

Returns
0 if success
i if the QR algorithm failed to compute all the eigenvalues 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,out]optsOptions.
  • Output parameters opts.n_aed, opts.n_sweep and opts.n_shifts_total are updated inside the routine.

◆ multishift_QR_sweep()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, bool > = true>
void tlapack::multishift_QR_sweep ( bool  want_t,
bool  want_z,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
const vector_t s,
matrix_t Z 
)

multishift_QR_sweep performs a single small-bulge multi-shift QR sweep.

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. ilo and ihi determine an isolated block in A.
[in,out]An by n matrix. Hessenberg matrix on which AED will be performed
[in]scomplex vector. Vector containing the shifts to be used during the sweep
[in,out]Zn by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A accumulated into Z.

◆ multishift_QZ_sweep()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR alpha_t, TLAPACK_VECTOR beta_t, enable_if_t< is_complex< type_t< alpha_t > >, bool > = true>
void tlapack::multishift_QZ_sweep ( bool  want_t,
bool  want_q,
bool  want_z,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
matrix_t A,
matrix_t B,
const alpha_t alpha,
const beta_t beta,
matrix_t Q,
matrix_t Z 
)

multishift_QR_sweep performs a single small-bulge multi-shift QR sweep.

Parameters
[in]want_tbool. If true, the full Schur factor T will be computed.
[in]want_qbool. If true, the Schur vectors Q 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. ilo and ihi determine an isolated block in (A,B).
[in]An by n matrix. Hessenberg matrix on which AED will be performed
[in]Bn by n matrix. Hessenberg matrix on which AED will be performed
[in]alphacomplex vector. Vector containing the shifts to be used during the sweep
[in]betavector. Vector containing the scale factor of the shifts to be used during the sweep
[in]Qn by n matrix. On entry, the previously calculated Schur factors.
[in]Zn by n matrix. On entry, the previously calculated Schur factors.

◆ ung2l()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ung2l ( matrix_t A,
const vector_t tau 
)

Generates an m-by-n matrix Q with orthonormal columns, which is defined as the last n columns of a product of k elementary reflectors of order m.

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

The reflectors are stored in the matrix A as returned by geqlf

Parameters
[in,out]Am-by-n matrix. On entry, the (n+k-i)-th column must contains the vector which defines the elementary reflector \(H_i\), for \(i=0,1,...,k-1\), as returned by geqlf. On exit, the m-by-n matrix \(Q\).
[in]tauReal vector of length min(m,n). The scalar factors of the elementary reflectors.
Returns
0 if success

◆ ung2r()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ung2r ( matrix_t A,
const vector_t tau 
)

Generates a matrix Q with orthogonal columns.

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

Parameters
[in,out]Am-by-n matrix. On entry, the i-th column must contains the vector which defines the elementary reflector \(H_i\), for \(i=0,1,...,k-1\), as returned by geqrf. On exit, the m-by-n matrix \(Q\).
[in]tauReal vector of length min(m,n). The scalar factors of the elementary reflectors.
Returns
0 if success

◆ ungbr_p()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::ungbr_p ( const size_type< matrix_t k,
matrix_t A,
const vector_t tau,
const UngbrOpts opts = {} 
)

Generates the unitary matrix P**H determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H.

P**H is defined as a product of elementary reflectors G(i).

A is assumed to have been an K-by-N matrix, and P**H is of order N: if k < n, P**H = G(k) . . . G(2) G(1) and ungbr_p returns the first m rows of P**H, where n >= m >= k; if k >= n, P**H = G(n-1) . . . G(2) G(1) and ungbr_p returns P**H as an N-by-N matrix.

Parameters
[in]kinteger. k is the number of rows in the original k-by-n matrix reduced by gebrd.
[in,out]Am-by-n matrix. On entry, the vectors which define the elementary reflectors as returned by gebrd. On exit, the m-by-n matrix P**H.
[in]tauvector. tau is a vector of length min(k,n) tau(i) must contain the scalar factor of the elementary reflector G(i), which determines P**H, as returned by gebrd in its array argument taup.
[in]optsOptions.

◆ ungbr_q()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::ungbr_q ( const size_type< matrix_t k,
matrix_t A,
const vector_t tau,
const UngbrOpts opts = {} 
)

Generates the unitary matrix Q determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H.

Q is defined as a product of elementary reflectors H(i).

A is assumed to have been an M-by-K matrix, and Q is of order M: if m >= k, Q = H(1) H(2) . . . H(k) and ungbr returns the first n columns of Q, where m >= n >= k; if m < k, Q = H(1) H(2) . . . H(m-1) and ungbr returns Q as an M-by-M matrix.

Parameters
[in]kinteger. k is the number of columns in the original m-by-k matrix reduced by gebrd.
[in,out]Am-by-n matrix. On entry, the vectors which define the elementary reflectors as returned by gebrd. On exit, the m-by-n matrix Q.
[in]tauvector. tau is a vector of length min(m,k) tau(i) must contain the scalar factor of the elementary reflector H(i), which determines Q, as returned by gebrd in its array argument tauq.
[in]optsOptions.

◆ ungl2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ungl2 ( matrix_t Q,
const vector_t tauw 
)

Generates all or part of the unitary matrix Q from an LQ factorization determined by gelq2 (unblocked algorithm).

The matrix Q is defined as the first k rows of a product of k elementary reflectors of order n

\[ Q = H(k)**H ... H(2)**H H(1)**H \]

as returned by gelq2 and k <= n.

Returns
0 if success
Parameters
[in,out]Qk-by-n matrix. On entry, the i-th row must contain the vector which defines the elementary reflector H(j), for j = 1,2,...,k, as returned by gelq2 in the first k rows of its array argument A. On exit, the k by n matrix Q.
[in]tauwComplex vector of length min(m,n). tauw(j) must contain the scalar factor of the elementary reflector H(j), as returned by gelq2.

◆ ungq()

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

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

Blocked algorithm.

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:
See also
ungq_level2_work().
Parameters
[in,out]Am-by-n matrix. On entry,
[in]tauVector of length k. Scalar factors of the elementary reflectors.
[in]optsOptions.
Returns
0 if success.

◆ ungq_level2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::ungq_level2 ( direction_t  direction,
storage_t  storeMode,
matrix_t A,
const vector_t tau 
)

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:
See also
ungq_level2_work().
Parameters
[in,out]Am-by-n matrix. On entry,
[in]tauVector of length k. Scalar factors of the elementary reflectors.
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 )

◆ ungr2()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t>
int tlapack::ungr2 ( matrix_t A,
const vector_t tau 
)

Generates an m-by-n matrix Q with orthonormal columns, which is defined as the last m rows of a product of k elementary reflectors of order n.

\[ Q = H_1^H H_2^H ... H_k^H \]

The reflectors are stored in the matrix A as returned by gerqf

Parameters
[in,out]Am-by-n matrix. On entry, the (m-k+i)-th row must contain the vector which defines the elementary reflector H(i), for i = 1,2,...,k, as returned by GERQF in the last k rows of its matrix argument A. On exit, the m-by-n matrix \(Q\).
[in]tauReal vector of length min(m,n). The scalar factors of the elementary reflectors.
Returns
0 if success

◆ unm2l()

template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unm2l ( side_t  side,
trans_t  trans,
const matrixA_t A,
const tau_t tau,
matrixC_t C 
)

Applies unitary matrix Q from an QL factorization to a matrix C.

The matrix Q is represented as a product of elementary reflectors as returned by geqlf

Template Parameters
side_tEither Side or any class that implements operator Side().
trans_tEither Op or any class that implements operator Op().
Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]A
[in]tauVector of length k Contains the scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:

◆ unm2r()

template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unm2r ( side_t  side,
trans_t  trans,
const matrixA_t A,
const tau_t tau,
matrixC_t C 
)

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].

Template Parameters
side_tEither Side or any class that implements operator Side().
trans_tEither Op or any class that implements operator Op().
Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]A
[in]tauVector of length k Contains the scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:

◆ unmhr()

template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t>
int tlapack::unmhr ( Side  side,
Op  trans,
size_type< matrix_t ilo,
size_type< matrix_t ihi,
const matrix_t A,
const vector_t tau,
matrix_t C 
)

Applies unitary matrix Q to a matrix C.

Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[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]An-by-n matrix Matrix containing orthogonal vectors, as returned by gehrd
[in]tauVector of length n-1 Contains the scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:

◆ unml2()

template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unml2 ( side_t  side,
trans_t  trans,
const matrixA_t A,
const tau_t tau,
matrixC_t C 
)

Applies unitary matrix Q from an LQ factorization to a matrix C.

The matrix Q is represented as a product of elementary reflectors as returned by gelqf

Template Parameters
side_tEither Side or any class that implements operator Side().
trans_tEither Op or any class that implements operator Op().
Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]A
[in]tauVector of length k Contains the scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following:

◆ unmq()

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::unmq ( side_t  side,
trans_t  trans,
direction_t  direction,
storage_t  storeMode,
const matrixV_t V,
const vector_t tau,
matrixC_t C,
const UnmqOpts opts = {} 
)

Applies unitary matrix Q to a matrix C.

Blocked algorithm.

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.
Returns
0 if success.

◆ unmq_level2()

template<TLAPACK_SMATRIX matrixV_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR vector_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t>
int tlapack::unmq_level2 ( side_t  side,
trans_t  trans,
direction_t  direction,
storage_t  storeMode,
const matrixV_t V,
const vector_t tau,
matrixC_t C 
)

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:
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 )

◆ unmr2()

template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_VECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t>
int tlapack::unmr2 ( side_t  side,
trans_t  trans,
const matrixA_t A,
const tau_t tau,
matrixC_t C 
)

Applies unitary matrix Q from an RQ factorization to a matrix C.

The matrix Q is represented as a product of elementary reflectors as returned by gerqf

Template Parameters
side_tEither Side or any class that implements operator Side().
trans_tEither Op or any class that implements operator Op().
Parameters
[in]sideSpecifies which side op(Q) is to be applied.
[in]transThe operation \(op(Q)\) to be used:
[in]A
[in]tauVector of length k Contains the scalar factors of the elementary reflectors.
[in,out]Cm-by-n matrix. On exit, C is replaced by one of the following: