<T>LAPACK 0.1.2
C++ Template Linear Algebra PACKage
|
All other computational routines that are not part of the BLAS. More...
Functions | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR alpha_t, TLAPACK_SVECTOR beta_t> | |
void | tlapack::aggressive_early_deflation_generalized (bool want_s, bool want_q, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, size_type< matrix_t > nw, matrix_t &A, matrix_t &B, alpha_t &alpha, beta_t &beta, matrix_t &Q, matrix_t &Z, size_type< matrix_t > &ns, size_type< matrix_t > &nd, FrancisOpts &opts) |
aggressive_early_deflation_generalized accepts as input an upper Hessenberg pencil (A,B) and performs a unitary similarity transformation designed to detect and deflate fully converged eigenvalues from a trailing principal subpencil. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t, enable_if_t< is_complex< type_t< vector_t > >, int > > | |
void | tlapack::aggressive_early_deflation_work (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, work_t &work, 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. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gebd2_work (matrix_t &A, vector_t &tauv, vector_t &tauw, work_t &work) |
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_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gebrd_work (matrix_t &A, vector_t &tauv, vector_t &tauw, work_t &work, const GebrdOpts &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_WORKSPACE work_t> | |
int | tlapack::gehd2_work (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau, work_t &work) |
Reduces a general square matrix to upper Hessenberg form. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gehrd_work (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, vector_t &tau, work_t &work, const GehrdOpts &opts={}) |
Reduces a general square matrix to upper Hessenberg form. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gelq2_work (matrix_t &A, vector_t &tauw, work_t &work) |
Computes an LQ factorization of a complex m-by-n matrix A using an unblocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gelqf_work (A_t &A, tau_t &tau, work_t &work, const GelqfOpts &opts={}) |
Computes an LQ factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gelqt_work (matrix_t &A, matrix_t &TT, work_t &work) |
Computes an LQ factorization of a complex m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::geql2_work (matrix_t &A, vector_t &tau, work_t &work) |
Computes a QL factorization of a matrix A. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::geqlf_work (A_t &A, tau_t &tau, work_t &work, const GeqlfOpts &opts={}) |
Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::geqr2_work (matrix_t &A, vector_t &tau, work_t &work) |
Computes a QR factorization of a matrix A. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::geqrf_work (A_t &A, tau_t &tau, work_t &work, const GeqrfOpts &opts={}) |
Computes a QR factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gerq2_work (matrix_t &A, vector_t &tau, work_t &work) |
Computes an RQ factorization of a matrix A. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::gerqf_work (A_t &A, tau_t &tau, work_t &work, const GerqfOpts &opts={}) |
Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR r_vector_t> | |
int | tlapack::gesvd (bool want_u, bool want_vt, matrix_t &A, r_vector_t &s, matrix_t &U, matrix_t &Vt, const GesvdOpts &opts={}) |
Computes the singular values and, optionally, the right and/or left singular vectors from the singular value decomposition (SVD) of a real M-by-N matrix A. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR piv_t> | |
int | tlapack::getrf_level0 (matrix_t &A, piv_t &piv) |
getrf computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR piv_t> | |
int | tlapack::getrf_recursive (matrix_t &A, piv_t &piv) |
getrf_recursive computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges. | |
template<TLAPACK_MATRIX matrix_t> | |
int | tlapack::getri_uili (matrix_t &A) |
getri_uili calculates the inverse of a general n-by-n matrix A | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::getri_uxli_work (matrix_t &A, work_t &work) |
getri computes inverse of a general n-by-n matrix A by solving for X in the following equation Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SMATRIX B_t, TLAPACK_SMATRIX Q_t, TLAPACK_SMATRIX Z_t> | |
int | tlapack::gghd3 (bool wantq, bool wantz, size_type< A_t > ilo, size_type< A_t > ihi, A_t &A, B_t &B, Q_t &Q, Z_t &Z, const Gghd3Opts &opts={}) |
Reduces a pair of real square matrices (A, B) to generalized upper Hessenberg form using unitary transformations, where A is a general matrix and B is upper triangular. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SMATRIX B_t, TLAPACK_SMATRIX Q_t, TLAPACK_SMATRIX Z_t> | |
int | tlapack::gghrd (bool wantq, bool wantz, size_type< A_t > ilo, size_type< A_t > ihi, A_t &A, B_t &B, Q_t &Q, Z_t &Z) |
Reduces a pair of real square matrices (A, B) to generalized upper Hessenberg form using unitary transformations, where A is a general matrix and B is upper triangular. | |
template<TLAPACK_SMATRIX T_t, TLAPACK_SVECTOR CL_t, TLAPACK_SVECTOR SL_t, TLAPACK_SVECTOR CR_t, TLAPACK_SVECTOR SR_t> | |
void | tlapack::hessenberg_rq (T_t &T, CL_t &cl, SL_t &sl, CR_t &cr, SR_t &sr) |
Applies a sequence of rotations to an upper triangular matrix T from the left (making it an upper Hessenberg matrix) and reduces that matrix back to upper triangular form using rotations from the right. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_SVECTOR tau_t, class uplo_t > | |
int | tlapack::hetd2 (uplo_t uplo, A_t &A, tau_t &tau) |
Reduce a hermitian matrix to real symmetric tridiagonal form by a unitary similarity transformation: Q**H * A * Q = T. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t, enable_if_t< is_complex< type_t< vector_t > >, bool > = true> | |
void | tlapack::multishift_QR_sweep_work (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, work_t &work) |
multishift_QR_sweep performs a single small-bulge multi-shift QR sweep. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t, enable_if_t< is_complex< type_t< vector_t > >, int > = 0> | |
int | tlapack::multishift_qr_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, 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. Workspace is provided as an argument. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t, disable_if_allow_optblas_t< matrix_t > = 0> | |
int | tlapack::potf2 (uplo_t uplo, matrix_t &A) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a level-2 algorithm. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
int | tlapack::potrf2 (uplo_t uplo, matrix_t &A, const EcOpts &opts={}) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using the recursive algorithm. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
int | tlapack::potrf_blocked (uplo_t uplo, matrix_t &A, const BlockedCholeskyOpts &opts) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a blocked algorithm. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
int | tlapack::potrf_rl (uplo_t uplo, matrix_t &A, const BlockedCholeskyOpts &opts={}) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a blocked algorithm. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrixA_t, TLAPACK_MATRIX matrixB_t> | |
int | tlapack::potrs (uplo_t uplo, const matrixA_t &A, matrixB_t &B) |
Apply the Cholesky factorization to solve a linear system. | |
template<TLAPACK_SIDE side_t, TLAPACK_DIRECTION direction_t, TLAPACK_SVECTOR C_t, TLAPACK_SVECTOR S_t, TLAPACK_SMATRIX A_t> | |
int | tlapack::rot_sequence (side_t side, direction_t direction, const C_t &c, const S_t &s, A_t &A) |
Applies a sequence of plane rotations to an (m-by-n) matrix. | |
template<TLAPACK_SIDE side_t, TLAPACK_DIRECTION direction_t, TLAPACK_SMATRIX C_t, TLAPACK_SMATRIX S_t, TLAPACK_SMATRIX A_t> | |
void | tlapack::rot_sequence3 (side_t side, direction_t direction, const C_t &C, const S_t &S, A_t &A) |
Applies a sequence of plane rotations to an (m-by-n) matrix. | |
template<class matrix_t , class d_t , class e_t , enable_if_t< is_same_v< type_t< d_t >, real_type< type_t< d_t > > >, int > = 0, enable_if_t< is_same_v< type_t< e_t >, real_type< type_t< e_t > > >, int > = 0> | |
int | tlapack::svd_qr (Uplo uplo, bool want_u, bool want_vt, d_t &d, e_t &e, matrix_t &U, matrix_t &Vt) |
Computes the singular values and, optionally, the right and/or left singular vectors from the singular value decomposition (SVD) of a real N-by-N (upper or lower) bidiagonal matrix B using the implicit zero-shift QR algorithm. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
int | tlapack::trtri_recursive (uplo_t uplo, Diag diag, matrix_t &C, const EcOpts &opts={}) |
TRTRI computes the inverse of a triangular matrix in-place Input is a triangular matrix, output is its inverse This is the recursive variant. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ung2l_work (matrix_t &A, const vector_t &tau, work_t &work) |
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, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ung2r_work (matrix_t &A, const vector_t &tau, work_t &work) |
Generates a matrix Q with orthogonal columns. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungbr_p_work (const size_type< matrix_t > k, matrix_t &A, const vector_t &tau, work_t &work, 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. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungbr_q_work (const size_type< matrix_t > k, matrix_t &A, const vector_t &tau, work_t &work, const UngbrOpts &opts={}) |
Generates the unitary matrix Q determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::unghr_work (size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, const vector_t &tau, work_t &work) |
Generates a m-by-n matrix Q with orthogonal columns. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungl2_work (matrix_t &Q, const vector_t &tauw, work_t &work) |
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> | |
int | tlapack::unglq (matrix_t &A, const vector_t &tau, const UnglqOpts &opts={}) |
Generates all or part of the unitary matrix Q from an LQ factorization determined by gelqf. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungq_level2_work (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau, work_t &work) |
Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungq_work (direction_t direction, storage_t storeMode, matrix_t &A, const vector_t &tau, work_t &work, const UngqOpts &opts={}) |
Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t> | |
int | tlapack::ungql (matrix_t &A, const vector_t &tau, const UngqlOpts &opts={}) |
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_SVECTOR vector_t> | |
int | tlapack::ungqr (matrix_t &A, const vector_t &tau, const UngqrOpts &opts={}) |
Generates a matrix Q with orthogonal columns. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::ungr2_work (matrix_t &A, const vector_t &tau, work_t &work) |
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 matrix_t, TLAPACK_SVECTOR vector_t> | |
int | tlapack::ungrq (matrix_t &A, const vector_t &tau, const UngrqOpts &opts={}) |
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 Q_t, TLAPACK_SVECTOR tau_t, class uplo_t > | |
int | tlapack::ungtr (uplo_t uplo, Q_t &Q, const tau_t &tau) |
Generates a real orthogonal matrix Q which is defined as the product of k elementary reflectors of order n, as returned by hetrd. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t> | |
int | tlapack::unmhr_work (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, work_t &work) |
Applies unitary matrix Q to a matrix C. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t> | |
int | tlapack::unmlq (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C, const UnmlqOpts &opts={}) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code. | |
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, TLAPACK_WORKSPACE work_t> | |
int | tlapack::unmq_level2_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) |
Applies unitary matrix Q to a matrix C. 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, TLAPACK_WORKSPACE work_t> | |
int | tlapack::unmq_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 UnmqOpts &opts={}) |
Applies unitary matrix Q to a matrix C. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t> | |
int | tlapack::unmql (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C, const UnmqlOpts &opts={}) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t> | |
int | tlapack::unmqr (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C, const UnmqrOpts &opts={}) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SVECTOR tau_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t> | |
int | tlapack::unmrq (side_t side, trans_t trans, const matrixA_t &A, const tau_t &tau, matrixC_t &C, const UnmrqOpts &opts={}) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code. | |
All other computational routines that are not part of the BLAS.
void tlapack::aggressive_early_deflation_generalized | ( | bool | want_s, |
bool | want_q, | ||
bool | want_z, | ||
size_type< matrix_t > | ilo, | ||
size_type< matrix_t > | ihi, | ||
size_type< matrix_t > | nw, | ||
matrix_t & | A, | ||
matrix_t & | B, | ||
alpha_t & | alpha, | ||
beta_t & | beta, | ||
matrix_t & | Q, | ||
matrix_t & | Z, | ||
size_type< matrix_t > & | ns, | ||
size_type< matrix_t > & | nd, | ||
FrancisOpts & | opts | ||
) |
aggressive_early_deflation_generalized accepts as input an upper Hessenberg pencil (A,B) and performs a unitary similarity transformation designed to detect and deflate fully converged eigenvalues from a trailing principal subpencil.
On output (A,B) has been over- written by a new perturbation that is a perturbation of an orthogonal similarity transformation of (A,B). It is to be hoped that the final version of (A,B) has many zero subdiagonal entries.
[in] | want_s | bool. If true, the full Schur factors (H,T) will be computed. |
[in] | want_q | bool. If true, the Schur vectors Q will be computed. |
[in] | want_z | bool. If true, the Schur vectors Z will be computed. |
[in] | ilo | integer. Either ilo=0 or A(ilo,ilo-1) = 0. |
[in] | ihi | integer. ilo and ihi determine an isolated block in A. |
[in] | nw | integer. 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] | A | n by n matrix. Hessenberg matrix on which AED will be performed |
[in,out] | B | n by n matrix. Upper triangular matrix on which AED will be performed |
[out] | alpha | size n vector. |
[out] | beta | size n vector. On exit, the entries (alpha,beta)[ihi-nd-ns:ihi-nd] contain the unconverged eigenvalues that can be used a shifts. The entries (alpha,beta)[ihi-nd:ihi] contain the converged eigenvalues. Entries outside the range [ihi-nw:ihi] are not changed. The converged shifts are stored in the same positions as their correspinding diagonal elements in A. |
[in,out] | Q | n by n matrix. On entry, the previously calculated Schur factors On exit, the left orthogonal updates applied to (A,B) accumulated into Q. |
[in,out] | Z | n by n matrix. On entry, the previously calculated Schur factors On exit, the right orthogonal updates applied to (A,B) accumulated into Z. |
[out] | ns | integer. Number of eigenvalues available as shifts in s. |
[out] | nd | integer. Number of converged eigenvalues available as shifts in s. |
[in,out] | opts | Options.
|
void tlapack::aggressive_early_deflation_work | ( | 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, | ||
work_t & | work, | ||
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. Workspace is provided as an argument.
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.
[in] | want_t | bool. If true, the full Schur factor T will be computed. |
[in] | want_z | bool. If true, the Schur vectors Z will be computed. |
[in] | ilo | integer. Either ilo=0 or A(ilo,ilo-1) = 0. |
[in] | ihi | integer. ilo and ihi determine an isolated block in A. |
[in] | nw | integer. 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] | A | n by n matrix. Hessenberg matrix on which AED will be performed |
[out] | s | size 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] | Z | n by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A accumulated into Z. |
[out] | ns | integer. Number of eigenvalues available as shifts in s. |
[out] | nd | integer. Number of converged eigenvalues available as shifts in s. |
[in,out] | opts | Options.
|
work | Workspace. Use the workspace query to determine the size needed. |
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 * 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).
[in,out] | A | m-by-n matrix. On entry, the m by n general matrix to be reduced.
|
[out] | tauv | Real vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Q. |
[out] | tauw | Real vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Z. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gebrd_work | ( | matrix_t & | A, |
vector_t & | tauv, | ||
vector_t & | tauw, | ||
work_t & | work, | ||
const GebrdOpts & | 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, \]
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).
[in,out] | A | m-by-n matrix. On entry, the m by n general matrix to be reduced. On exit,
|
[out] | tauv | vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix Q. |
[out] | tauw | vector of length min(m,n). The scalar factors of the elementary reflectors which represent the unitary matrix P. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gehd2_work | ( | size_type< matrix_t > | ilo, |
size_type< matrix_t > | ihi, | ||
matrix_t & | A, | ||
vector_t & | tau, | ||
work_t & | work | ||
) |
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].
[in] | ilo | integer |
[in] | ihi | integer 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] | A | n-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] | tau | Real vector of length n-1. The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gehrd_work | ( | size_type< matrix_t > | ilo, |
size_type< matrix_t > | ihi, | ||
matrix_t & | A, | ||
vector_t & | tau, | ||
work_t & | work, | ||
const GehrdOpts & | 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].
[in] | ilo | integer |
[in] | ihi | integer 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] | A | n-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] | tau | Real vector of length n-1. The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes an LQ factorization of a complex m-by-n matrix A using an unblocked algorithm. Workspace is provided as an argument.
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].
[in,out] | A | m-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] | tauw | Complex vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes an LQ factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument.
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.
[in,out] | A | m-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] | tau | min(n,m) vector. The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gelqt_work | ( | matrix_t & | A, |
matrix_t & | TT, | ||
work_t & | work | ||
) |
Computes an LQ factorization of a complex m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument.
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).
[in,out] | A | m-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] | TT | m-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 |
work | Workspace. Use the workspace query to determine the size needed. |
Computes a QL factorization of a 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[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].
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm. Workspace is provided as an argument.
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.
[in,out] | A | m-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] | tau | min(n,m) vector. The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes a QR factorization of a 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].
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes a QR factorization of an m-by-n matrix A using a blocked algorithm. 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].
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes an RQ factorization of a 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[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].
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
Computes an RQ factorization of an m-by-n matrix A using a blocked algorithm. 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[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].
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gesvd | ( | bool | want_u, |
bool | want_vt, | ||
matrix_t & | A, | ||
r_vector_t & | s, | ||
matrix_t & | U, | ||
matrix_t & | Vt, | ||
const GesvdOpts & | opts = {} |
||
) |
Computes the singular values and, optionally, the right and/or left singular vectors from the singular value decomposition (SVD) of a real M-by-N matrix A.
The SVD of A has the form B = U * S * V^H where S is the diagonal matrix of singular values, U is a unitary matrix of left singular vectors, and V is a unitary matrix of right singular vectors. Depending on the dimensions of U and Vt, either the reduced or full unitary factors are determined.
[in] | want_u | bool |
[in] | want_vt | bool |
[in,out] | A | m-by-n matrix. |
[out] | s | vector of length min(m,n). The singular values of A, sorted so that S(i) >= S(i+1). |
[in,out] | U | m-by-m matrix. |
[in,out] | Vt | n-by-n matrix. |
[in] | opts | Options.
|
int tlapack::getrf_level0 | ( | matrix_t & | A, |
piv_t & | piv | ||
) |
getrf computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.
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).
This is a Level 0 version of the algorithm.
[in,out] | A | m-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] | piv | is 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 |
int tlapack::getrf_recursive | ( | matrix_t & | A, |
piv_t & | piv | ||
) |
getrf_recursive computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.
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).
This is a recursive version of the algorithm.
[in,out] | A | m-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] | piv | is 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 |
int tlapack::getri_uili | ( | matrix_t & | A | ) |
getri_uili calculates the inverse of a general n-by-n matrix A
A is assumed to be in the form L U factors on the input, trtri is used to invert U and L in place thereafter, ul_mult is called to calculate U^(-1)L^(-1) in place.
[in,out] | A | n-by-n matrix. On entry, the factors L and U from the factorization 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. |
int tlapack::getri_uxli_work | ( | matrix_t & | A, |
work_t & | work | ||
) |
getri computes inverse of a general n-by-n matrix A by solving for X in the following equation Workspace is provided as an argument.
\[ U X L = I \]
[in,out] | A | n-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. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::gghd3 | ( | bool | wantq, |
bool | wantz, | ||
size_type< A_t > | ilo, | ||
size_type< A_t > | ihi, | ||
A_t & | A, | ||
B_t & | B, | ||
Q_t & | Q, | ||
Z_t & | Z, | ||
const Gghd3Opts & | opts = {} |
||
) |
Reduces a pair of real square matrices (A, B) to generalized upper Hessenberg form using unitary transformations, where A is a general matrix and B is upper triangular.
[in] | wantq | boolean |
[in] | wantz | boolean |
[in] | ilo | integer |
[in] | ihi | integer |
[in,out] | A | n-by-n matrix. |
[in,out] | B | n-by-n matrix. |
[in,out] | Q | n-by-n matrix. |
[in,out] | Z | n-by-n matrix. |
[in] | opts | Options. |
int tlapack::gghrd | ( | bool | wantq, |
bool | wantz, | ||
size_type< A_t > | ilo, | ||
size_type< A_t > | ihi, | ||
A_t & | A, | ||
B_t & | B, | ||
Q_t & | Q, | ||
Z_t & | Z | ||
) |
Reduces a pair of real square matrices (A, B) to generalized upper Hessenberg form using unitary transformations, where A is a general matrix and B is upper triangular.
[in] | wantq | boolean |
[in] | wantz | boolean |
[in] | ilo | integer |
[in] | ihi | integer |
[in,out] | A | n-by-n matrix. |
[in,out] | B | n-by-n matrix. |
[in,out] | Q | n-by-n matrix. |
[in,out] | Z | n-by-n matrix. |
|
inline |
Applies a sequence of rotations to an upper triangular matrix T from the left (making it an upper Hessenberg matrix) and reduces that matrix back to upper triangular form using rotations from the right.
i.e. rot_sequence(LEFT_SIDE, FORWARD, cl, sl, T1) = rot_sequence(RIGHT_SIDE, FORWARD, cr, sr, T2)
This is an important component of gghd3
[in,out] | T | n-by-n Upper triangular matrix. |
[in] | cl | Real vector of length n-1. Cosines of the left rotations |
[in] | sl | Vector of length n-1. Sines of the left rotations |
[out] | cr | Real vector of length n-1. Cosines of the right rotations |
[out] | sr | Vector of length n-1. Sines of the right rotations |
Reduce a hermitian matrix to real symmetric tridiagonal form by a unitary similarity transformation: Q**H * A * Q = T.
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | A | n-by-n symmetric matrix. On exit, the main diagonal and offdiagonal contain the elements of the symmetric tridiagonal matrix B. The other positions are used to store elementary Householder reflectors. |
[out] | tau | Vector of length n-1. The scalar factors of the elementary reflectors. |
void tlapack::multishift_QR_sweep_work | ( | 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, | ||
work_t & | work | ||
) |
multishift_QR_sweep performs a single small-bulge multi-shift QR sweep. Workspace is provided as an argument.
[in] | want_t | bool. If true, the full Schur factor T will be computed. |
[in] | want_z | bool. If true, the Schur vectors Z will be computed. |
[in] | ilo | integer. Either ilo=0 or A(ilo,ilo-1) = 0. |
[in] | ihi | integer. ilo and ihi determine an isolated block in A. |
[in,out] | A | n by n matrix. Hessenberg matrix on which AED will be performed |
[in] | s | complex vector. Vector containing the shifts to be used during the sweep |
[in,out] | Z | n by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A accumulated into Z. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::multishift_qr_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, | ||
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. Workspace is provided as an argument.
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.
[in] | want_t | bool. If true, the full Schur factor T will be computed. |
[in] | want_z | bool. If true, the Schur vectors Z will be computed. |
[in] | ilo | integer. Either ilo=0 or A(ilo,ilo-1) = 0. |
[in] | ihi | integer. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n. |
[in,out] | A | n 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] | w | size 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] | Z | n by n matrix. On entry, the previously calculated Schur factors On exit, the orthogonal updates applied to A are accumulated into Z. |
[in,out] | opts | Options.
|
work | Workspace. Use the workspace query to determine the size needed. |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a level-2 algorithm.
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.
uplo_t | Access type: Upper or Lower. Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | A | On entry, the Hermitian matrix A. |
int tlapack::potrf2 | ( | uplo_t | uplo, |
matrix_t & | A, | ||
const EcOpts & | opts = {} |
||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using the recursive algorithm.
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.
This is the recursive version of the algorithm. It divides the matrix into four submatrices:
\[ A = \begin{bmatrix} A_{11} & A_{12} \\ A_{21} & A_{22} \end{bmatrix} \]
where \(A_{11}\) is n1-by-n1 and \(A_{22}\) is n2-by-n2, with n1 = n/2 and n2 = n-n1, where n is the order of the matrix A. The subroutine calls itself to factor \(A_{11},\) updates and scales \(A_{21}\) or \(A_{12},\) updates \(A_{22},\) and calls itself to factor \(A_{22}.\)
uplo_t | Access type: Upper or Lower. Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | A | On entry, the Hermitian matrix A. |
[in] | opts | Options. Define the behavior of Exception Handling. |
int tlapack::potrf_blocked | ( | uplo_t | uplo, |
matrix_t & | A, | ||
const BlockedCholeskyOpts & | opts | ||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a blocked algorithm.
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.
uplo_t | Access type: Upper or Lower. Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | A | On entry, the Hermitian matrix A of size n-by-n. |
[in] | opts | Options. Define the behavior of checks for NaNs. |
int tlapack::potrf_rl | ( | uplo_t | uplo, |
matrix_t & | A, | ||
const BlockedCholeskyOpts & | opts = {} |
||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a blocked algorithm.
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.
uplo_t | Access type: Upper or Lower. Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | A | On entry, the Hermitian matrix A of size n-by-n. |
[in] | opts | Options. Define the behavior of checks for NaNs. |
Apply the Cholesky factorization to solve a linear system.
\[ A X = B, \]
where \(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.
uplo_t | Access type: Upper or Lower. Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in] | A | The factor U or L from the Cholesky factorization of A. |
[in,out] | B | On entry, the matrix B. On exit, the matrix X. |
int tlapack::rot_sequence | ( | side_t | side, |
direction_t | direction, | ||
const C_t & | c, | ||
const S_t & | s, | ||
A_t & | A | ||
) |
Applies a sequence of plane rotations to an (m-by-n) matrix.
When side = Side::Left, the transformation takes the form
A := P*A
and when side = Side::Right, the transformation takes the form
A := A*P**H
where P is an orthogonal matrix consisting of a sequence of k plane rotations, with k = m-1 when side = Side::Left and k = n-1 when side = Side::Right.
When direction = Direction::Forward, then
P = P(k-1) * P(k-2) * ... * P(1) * P(0)
and when direction = Direction::Backward, then
P = P(0) * P(1) * ... * P(k-2) * P(k-1),
when P(i) is a plane rotation matrix defined as
P(i) = ( 1 ) ( ... ) ( ... ) ( c(i) s(i) ) ( -conj(s(i)) c(i) ) ( 1 ) ( ... ) ( 1 )
which only modifies rows/columns i and i + 1.
[in] | side | Specifies whether the plane rotation matrix P is applied to A on the left or the right
|
[in] | direction | Specifies whether P is a forward or backward sequence of plane rotations. |
[in] | c | Real vector of length k. Cosines of the rotations |
[in] | s | Vector of length k. Sines of the rotations |
[in,out] | A | m-by-n matrix. |
void tlapack::rot_sequence3 | ( | side_t | side, |
direction_t | direction, | ||
const C_t & | C, | ||
const S_t & | S, | ||
A_t & | A | ||
) |
Applies a sequence of plane rotations to an (m-by-n) matrix.
When side = Side::Left, the transformation takes the form
A := P*A
and when side = Side::Right, the transformation takes the form
A := A*P**H
where P is an orthogonal matrix consisting of a sequence of k*l plane rotations, with k = m-1 when side = Side::Left and k = n-1 when side = Side::Right and l is the number of columns of C and S.
When direction = Direction::Forward, then
P = P(k-1,l-1) * ... * P(0,l-1) * ... * P(k-1,l-2) * ... * P(0,l-2)
and when direction = Direction::Backward, then
P = P(0,l-1) * ... * P(k-1,l-1) * ... * P(0,l-2) * ... * P(k-1,l-2)
where P(i,j) is a plane rotation matrix defined as
P(i,j) = ( 1 ) ( ... ) ( ... ) ( C(i,j) S(i,j) ) ( -conj(S(i,j)) C(i,j) ) ( 1 ) ( ... ) ( 1 )
which only modifies rows/columns i and i + 1.
This function is a variant of rot_sequence, which also applies sequence of plane rotations to a matrix. This variant can been seen as multiple applications of rot_sequence, but more cache-efficient.
Note: One of the implicit blocking parameters in this routine is l. The routine will work blocks in A of size 2*l-by-nb, where nb is a blocking parameter. If l is large, it may be better to call this routine multiple times.
Reference: "Restructuring the Tridiagonal and Bidiagonal QR algorithms for Performance" F. G. Van Zee, R. A. Van de Geijn, G. Quintana-Orti
[in] | side | Specifies whether the plane rotation matrix P is applied to A on the left or the right
|
[in] | direction | Specifies whether P is a forward or backward sequence of plane rotations. |
[in] | C | Real k-by-l matrix. Cosines of the rotations |
[in] | S | k-by-l matrix. Sines of the rotations |
[in,out] | A | m-by-n matrix. Matrix that the rotations will be applied to. |
int tlapack::svd_qr | ( | Uplo | uplo, |
bool | want_u, | ||
bool | want_vt, | ||
d_t & | d, | ||
e_t & | e, | ||
matrix_t & | U, | ||
matrix_t & | Vt | ||
) |
Computes the singular values and, optionally, the right and/or left singular vectors from the singular value decomposition (SVD) of a real N-by-N (upper or lower) bidiagonal matrix B using the implicit zero-shift QR algorithm.
The SVD of B has the form B = Q * S * P**T where S is the diagonal matrix of singular values, Q is an orthogonal matrix of left singular vectors, and P is an orthogonal matrix of right singular vectors. If left singular vectors are requested, this subroutine actually returns U*Q instead of Q, and, if right singular vectors are requested, this subroutine returns P**T*VT instead of P**T, for given real input matrices U and VT. When U and VT are the orthogonal matrices that reduce a general matrix A to bidiagonal form: A = U*B*VT, as computed by gebrd, then A = (U*Q) * S * (P**T*VT) is the SVD of A.
See "Computing Small Singular Values of Bidiagonal Matrices With Guaranteed High Relative Accuracy," by J. Demmel and W. Kahan, LAPACK Working Note #3 (or SIAM J. Sci. Statist. Comput. vol. 11, no. 5, pp. 873-912, Sept 1990) and "Accurate singular values and differential qd algorithms," by B. Parlett and V. Fernando, Technical Report CPAM-554, Mathematics Department, University of California at Berkeley, July 1992 for a detailed description of the algorithm.
[in] | uplo | Uplo::Upper, B is upper bidiagonal Uplo::Lower, B is lower bidiagonal |
[in] | want_u | bool |
[in] | want_vt | bool |
[in,out] | d | Real vector of length n. On entry, diagonal elements of the bidiagonal matrix B. On exit, the singular values of B in decreasing order. |
[in,out] | e | Real vector of length n-1. On entry, off-diagonal elements of the bidiagonal matrix B. On exit, the singular values of B in decreasing order. |
[in,out] | U | nu-by-m matrix. On entry, an nu-by-n unitary matrix. On exit, U is overwritten by U * Q. |
[in,out] | Vt | n-by-nvt matrix. On entry, an n-by-nvt unitary matrix. On exit, Vt is overwritten by P^H * Vt. |
int tlapack::trtri_recursive | ( | uplo_t | uplo, |
Diag | diag, | ||
matrix_t & | C, | ||
const EcOpts & | opts = {} |
||
) |
TRTRI computes the inverse of a triangular matrix in-place Input is a triangular matrix, output is its inverse This is the recursive variant.
[in] | uplo |
|
[in] | diag | Whether C has a unit or non-unit diagonal:
|
[in,out] | C | n-by-n matrix. On entry, the n-by-n triangular matrix to be inverted. On exit, the inverse. |
[in] | opts | Options. |
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
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
Generates a matrix Q with orthogonal columns.
\[ Q = H_1 H_2 ... H_k \]
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::ungbr_p_work | ( | const size_type< matrix_t > | k, |
matrix_t & | A, | ||
const vector_t & | tau, | ||
work_t & | work, | ||
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. Workspace is provided as an argument.
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.
[in] | k | integer. k is the number of rows in the original k-by-n matrix reduced by gebrd. |
[in,out] | A | m-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] | tau | vector. 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] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::ungbr_q_work | ( | const size_type< matrix_t > | k, |
matrix_t & | A, | ||
const vector_t & | tau, | ||
work_t & | work, | ||
const UngbrOpts & | opts = {} |
||
) |
Generates the unitary matrix Q determined by gebrd when reducing a matrix A to bidiagonal form: A = Q * B * P**H. Workspace is provided as an argument.
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.
[in] | k | integer. k is the number of columns in the original m-by-k matrix reduced by gebrd. |
[in,out] | A | m-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] | tau | vector. 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] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::unghr_work | ( | size_type< matrix_t > | ilo, |
size_type< matrix_t > | ihi, | ||
matrix_t & | A, | ||
const vector_t & | tau, | ||
work_t & | work | ||
) |
Generates a m-by-n matrix Q with orthogonal columns. Workspace is provided as an argument.
[in] | ilo | integer |
[in] | ihi | integer 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] | A | m-by-n matrix. On entry, the vectors which define the elementary reflectors. On exit, the m-by-n matrix Q. |
[in] | tau | Real vector of length n-1. The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
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.
[in,out] | Q | k-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] | tauw | Complex vector of length min(m,n). tauw(j) must contain the scalar factor of the elementary reflector H(j), as returned by gelq2. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::unglq | ( | matrix_t & | A, |
const vector_t & | tau, | ||
const UnglqOpts & | opts = {} |
||
) |
Generates all or part of the unitary matrix Q from an LQ factorization determined by gelqf.
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 gelqf and k <= n.
[in,out] | A | k-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] | tau | Complex vector of length min(m,n). tau(j) must contain the scalar factor of the elementary reflector H(j), as returned by gelqf. |
[in] | opts | Options. |
int tlapack::ungq_level2_work | ( | direction_t | direction, |
storage_t | storeMode, | ||
matrix_t & | A, | ||
const vector_t & | tau, | ||
work_t & | work | ||
) |
Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument.
[in] | direction | Indicates how Q is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in,out] | A | m-by-n matrix. On entry,
|
[in] | tau | Vector of length k. Scalar factors of the elementary reflectors. |
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 )
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::ungq_work | ( | direction_t | direction, |
storage_t | storeMode, | ||
matrix_t & | A, | ||
const vector_t & | tau, | ||
work_t & | work, | ||
const UngqOpts & | opts = {} |
||
) |
Generates a matrix Q that is the product of elementary reflectors. Workspace is provided as an argument.
Blocked algorithm.
[in] | direction | Indicates how Q is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in,out] | A | m-by-n matrix. On entry,
|
[in] | tau | Vector of length k. Scalar factors of the elementary reflectors. |
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::ungql | ( | matrix_t & | A, |
const vector_t & | tau, | ||
const UngqlOpts & | opts = {} |
||
) |
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
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
int tlapack::ungqr | ( | matrix_t & | A, |
const vector_t & | tau, | ||
const UngqrOpts & | opts = {} |
||
) |
Generates a matrix Q with orthogonal columns.
\[ Q = H_1 H_2 ... H_k \]
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
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
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::ungrq | ( | matrix_t & | A, |
const vector_t & | tau, | ||
const UngrqOpts & | opts = {} |
||
) |
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
[in,out] | A | m-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] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
[in] | opts | Options. |
Generates a real orthogonal matrix Q which is defined as the product of k elementary reflectors of order n, as returned by hetrd.
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in,out] | Q | n-by-n matrix. On entry, the vectors which define the elementary reflectors, as returned by hetrd. On exit, the n-by-n orthogonal matrix Q. |
[in] | tau | Vector of length k. The scalar factors of the elementary reflectors. |
int tlapack::unmhr_work | ( | 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, | ||
work_t & | work | ||
) |
Applies unitary matrix Q to a matrix C. Workspace is provided as an argument.
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | ilo | integer |
[in] | ihi | integer 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] | A | n-by-n matrix Matrix containing orthogonal vectors, as returned by gehrd |
[in] | tau | Vector of length n-1 Contains the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::unmlq | ( | side_t | side, |
trans_t | trans, | ||
const matrixA_t & | A, | ||
const tau_t & | tau, | ||
matrixC_t & | C, | ||
const UnmlqOpts & | opts = {} |
||
) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code.
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].
side_t | Either Side or any class that implements operator Side() . |
trans_t | Either Op or any class that implements operator Op() . |
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | A |
|
[in] | tau | Vector of length k Contains the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
[in] | opts | Options. |
int tlapack::unmq_level2_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 | ||
) |
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].
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | direction | Indicates how Q is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored:
|
[in] | V |
|
[in] | tau | Vector of length k. Scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
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 )
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::unmq_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 UnmqOpts & | opts = {} |
||
) |
Applies unitary matrix Q to a matrix C. Workspace is provided as an argument.
Blocked algorithm.
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | direction | Indicates how Q is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored:
|
[in] | V |
|
[in] | tau | Vector of length k. Scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
[in] | opts | Options. |
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::unmql | ( | side_t | side, |
trans_t | trans, | ||
const matrixA_t & | A, | ||
const tau_t & | tau, | ||
matrixC_t & | C, | ||
const UnmqlOpts & | opts = {} |
||
) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code.
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].
side_t | Either Side or any class that implements operator Side() . |
trans_t | Either Op or any class that implements operator Op() . |
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | A |
|
[in] | tau | Vector of length k Contains the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
[in] | opts | Options. |
int tlapack::unmqr | ( | side_t | side, |
trans_t | trans, | ||
const matrixA_t & | A, | ||
const tau_t & | tau, | ||
matrixC_t & | C, | ||
const UnmqrOpts & | opts = {} |
||
) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code.
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].
side_t | Either Side or any class that implements operator Side() . |
trans_t | Either Op or any class that implements operator Op() . |
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | A |
|
[in] | tau | Vector of length k Contains the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
[in] | opts | Options. |
int tlapack::unmrq | ( | side_t | side, |
trans_t | trans, | ||
const matrixA_t & | A, | ||
const tau_t & | tau, | ||
matrixC_t & | C, | ||
const UnmrqOpts & | opts = {} |
||
) |
Applies orthogonal matrix op(Q) to a matrix C using a blocked code.
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].
side_t | Either Side or any class that implements operator Side() . |
trans_t | Either Op or any class that implements operator Op() . |
[in] | side | Specifies which side op(Q) is to be applied.
|
[in] | trans | The operation \(op(Q)\) to be used:
|
[in] | A |
|
[in] | tau | Vector of length k Contains the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following:
|
[in] | opts | Options. |