<T>LAPACK 0.1.2
C++ Template Linear Algebra PACKage
|
Auxiliary routines that are used by the computational routines. More...
Classes | |
struct | tlapack::TestUploMatrix< T, idx_t, uplo, L > |
TestUploMatrix class. More... | |
Functions | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_generalized_similarity_transform (matrix_t &A, matrix_t &Q, matrix_t &Z, matrix_t &B) |
Calculates ||Q'*A*Z - B||. | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_generalized_similarity_transform (matrix_t &A, matrix_t &Q, matrix_t &Z, matrix_t &B, matrix_t &res, matrix_t &work) |
Calculates res = Q'*A*Z - B and the frobenius norm of res. | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_orthogonality (matrix_t &Q) |
Calculates ||Q'*Q - I||_F if m <= n or ||Q*Q' - I||_F otherwise. | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_orthogonality (matrix_t &Q, matrix_t &res) |
Calculates res = Q'*Q - I if m <= n or res = Q*Q' otherwise Also computes the frobenius norm of res. | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_similarity_transform (matrix_t &A, matrix_t &Q, matrix_t &B) |
Calculates ||Q'*A*Q - B||. | |
template<TLAPACK_MATRIX matrix_t> | |
real_type< type_t< matrix_t > > | tlapack::check_similarity_transform (matrix_t &A, matrix_t &Q, matrix_t &B, matrix_t &res, matrix_t &work) |
Calculates res = Q'*A*Q - B and the frobenius norm of res. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixB_t> | |
void | tlapack::conjtranspose (matrixA_t &A, matrixB_t &B, const TransposeOpts &opts={}) |
conjugate transpose a matrix A into a matrix B. | |
template<TLAPACK_VECTOR vector_t> | |
void | tlapack::conjugate (vector_t &x) |
Conjugates a vector. | |
template<TLAPACK_MATRIX matrix_t> | |
int | tlapack::generalized_schur_move (bool want_q, bool want_z, matrix_t &A, matrix_t &B, matrix_t &Q, matrix_t &Z, size_type< matrix_t > &ifst, size_type< matrix_t > &ilst) |
generalized_schur_move reorders the generalized Schur factorization of a pencil ( S, T ) = Q(A,B)Z**H so that the diagonal elements of (S,T) with row index IFST are moved to row ILST. | |
template<TLAPACK_CSMATRIX matrix_t, enable_if_t< is_real< type_t< matrix_t > >, bool > = true> | |
int | tlapack::generalized_schur_swap (bool want_q, bool want_z, matrix_t &A, matrix_t &B, matrix_t &Q, matrix_t &Z, const size_type< matrix_t > &j0, const size_type< matrix_t > &n1, const size_type< matrix_t > &n2) |
schur_swap, swaps 2 eigenvalues of A. | |
template<TLAPACK_MATRIX matrix_t> | |
auto | tlapack::infnorm_colmajor (const matrix_t &A) |
Calculates the infinity norm of a column-major matrix. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
auto | tlapack::infnorm_colmajor_work (const matrix_t &A, work_t &work) |
Calculates the infinity norm of a column-major matrix. Workspace is provided as an argument. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t> | |
auto | tlapack::infnorm_hermitian_colmajor (uplo_t uplo, const matrix_t &A) |
Calculates the infinity norm of a column-major hermitian matrix. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
auto | tlapack::infnorm_hermitian_colmajor_work (uplo_t uplo, const matrix_t &A, work_t &work) |
Calculates the infinity norm of a column-major hermitian matrix. Workspace is provided as an argument. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t> | |
auto | tlapack::infnorm_symmetric_colmajor (uplo_t uplo, const matrix_t &A) |
Calculates the infinity norm of a column-major symmetric matrix. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
auto | tlapack::infnorm_symmetric_colmajor_work (uplo_t uplo, const matrix_t &A, work_t &work) |
Calculates the infinity norm of a column-major symmetric matrix. Workspace is provided as an argument. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_DIAG diag_t, TLAPACK_MATRIX matrix_t> | |
auto | tlapack::infnorm_triangular_colmajor (uplo_t uplo, diag_t diag, const matrix_t &A) |
Calculates the infinity norm of a column-major triangular matrix. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_DIAG diag_t, TLAPACK_MATRIX matrix_t, TLAPACK_WORKSPACE work_t> | |
auto | tlapack::infnorm_triangular_colmajor_work (uplo_t uplo, diag_t diag, const matrix_t &A, work_t &work) |
Calculates the infinity norm of a column-major triangular matrix. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t> | |
void | tlapack::inv_house3 (const matrix_t &A, vector_t &v, type_t< vector_t > &tau) |
Inv_house calculates a reflector to reduce the first column in a 2x3 matrix A from the right to zero. | |
template<TLAPACK_SMATRIX A_t, TLAPACK_VECTOR vector_t, TLAPACK_SMATRIX X_t, TLAPACK_SMATRIX matrixY_t> | |
int | tlapack::labrd (A_t &A, vector_t &tauq, vector_t &taup, X_t &X, matrixY_t &Y) |
Reduces the first nb rows and columns of a general m by n matrix A to upper or lower bidiagonal form by an unitary transformation Q**H * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrixA_t, TLAPACK_MATRIX matrixB_t> | |
void | tlapack::lacpy (uplo_t uplo, const matrixA_t &A, matrixB_t &B) |
Copies a matrix from A to B. | |
template<TLAPACK_REAL real_t, enable_if_t<(is_real< real_t >), int > = 0> | |
void | tlapack::ladiv (const real_t &a, const real_t &b, const real_t &c, const real_t &d, real_t &p, real_t &q) |
Performs complex division in real arithmetic. | |
template<TLAPACK_COMPLEX T, enable_if_t< is_complex< T >, int > = 0> | |
T | tlapack::ladiv (const T &x, const T &y) |
Performs complex division in real arithmetic with complex arguments. | |
template<TLAPACK_CSMATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< vector_t > >, bool > = true, enable_if_t< is_real< type_t< matrix_t > >, bool > = true> | |
int | tlapack::lahqr (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) |
lahqr computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix, using the double-shift implicit QR algorithm. | |
template<TLAPACK_SCALAR T> | |
void | tlapack::lahqr_eig22 (T a00, T a01, T a10, T a11, complex_type< T > &s1, complex_type< T > &s2) |
Computes the eigenvalues of a 2x2 matrix A. | |
template<TLAPACK_REAL T> | |
void | tlapack::lahqr_schur22 (T &a, T &b, T &c, T &d, complex_type< T > &s1, complex_type< T > &s2, T &cs, T &sn) |
Computes the Schur factorization of a 2x2 matrix A. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_real< type_t< matrix_t > >, bool > = true> | |
int | tlapack::lahqr_shiftcolumn (const matrix_t &H, vector_t &v, complex_type< type_t< matrix_t > > s1, complex_type< type_t< matrix_t > > s2) |
Given a 2-by-2 or 3-by-3 matrix H, lahqr_shiftcolumn calculates a multiple of the product: (H - s1*I)*(H - s2*I)*e1. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, enable_if_t< is_complex< type_t< matrix_t > >, bool > = true> | |
int | tlapack::lahqr_shiftcolumn (const matrix_t &H, vector_t &v, type_t< matrix_t > s1, type_t< matrix_t > s2) |
Given a 2-by-2 or 3-by-3 matrix H, lahqr_shiftcolumn calculates a multiple of the product: (H - s1*I)*(H - s2*I)*e1. | |
template<TLAPACK_CSMATRIX matrix_t, TLAPACK_VECTOR alpha_t, TLAPACK_VECTOR beta_t> | |
int | tlapack::lahqz (bool want_s, bool want_q, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, matrix_t &B, alpha_t &alpha, beta_t &beta, matrix_t &Q, matrix_t &Z) |
lahqz computes the eigenvalues of a matrix pair (H,T), where H is an upper Hessenberg matrix and T is upper triangular, using the single/double-shift implicit QZ algorithm. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_VECTOR vector_t, bool > | |
int | tlapack::lahqz_shiftcolumn (const matrix_t &A, const matrix_t &B, vector_t &v, complex_type< type_t< matrix_t > > s1, complex_type< type_t< matrix_t > > s2, type_t< matrix_t > beta1, type_t< matrix_t > beta2) |
Given a 2-by-2 or 3-by-3 matrix pencil (A,B), lahqz_shiftcolumn calculates a multiple of the product: (beta2*A - s2*B)*B^(-1)*(beta1*A - s1*B)*B^(-1)*e1. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_VECTOR vector_t, TLAPACK_SMATRIX matrixT_t, TLAPACK_SMATRIX matrixY_t> | |
int | tlapack::lahr2 (size_type< matrix_t > k, size_type< matrix_t > nb, matrix_t &A, vector_t &tau, matrixT_t &T, matrixY_t &Y) |
Reduces a general square matrix to upper Hessenberg form. | |
template<TLAPACK_NORM norm_t, TLAPACK_SMATRIX matrix_t> | |
auto | tlapack::lange (norm_t normType, const matrix_t &A) |
Calculates the norm of a matrix. | |
template<TLAPACK_NORM norm_t, TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
auto | tlapack::lanhe (norm_t normType, uplo_t uplo, const matrix_t &A) |
Calculates the norm of a hermitian matrix. | |
template<TLAPACK_NORM norm_t, TLAPACK_UPLO uplo_t, TLAPACK_SMATRIX matrix_t> | |
auto | tlapack::lansy (norm_t normType, uplo_t uplo, const matrix_t &A) |
Calculates the norm of a symmetric matrix. | |
template<TLAPACK_NORM norm_t, TLAPACK_UPLO uplo_t, TLAPACK_DIAG diag_t, TLAPACK_SMATRIX matrix_t> | |
auto | tlapack::lantr (norm_t normType, uplo_t uplo, diag_t diag, const matrix_t &A) |
Calculates the norm of a symmetric matrix. | |
template<TLAPACK_REAL TX, TLAPACK_REAL TY, enable_if_t<(is_real< TX > &&is_real< TY >), int > = 0> | |
real_type< TX, TY > | tlapack::lapy2 (const TX &x, const TY &y) |
Finds \(\sqrt{x^2+y^2}\), taking care not to cause unnecessary overflow. | |
template<TLAPACK_REAL TX, TLAPACK_REAL TY, TLAPACK_REAL TZ, enable_if_t<(is_real< TX > &&is_real< TY > &&is_real< TZ >), int > = 0> | |
real_type< TX, TY, TZ > | tlapack::lapy3 (const TX &x, const TY &y, const TZ &z) |
Finds \(\sqrt{x^2+y^2+z^2}\), taking care not to cause unnecessary overflow or unnecessary underflow. | |
template<TLAPACK_SIDE side_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_SVECTOR vector_t, TLAPACK_SCALAR tau_t, TLAPACK_SMATRIX matrix_t, enable_if_t< std::is_convertible_v< direction_t, Direction >, int > = 0> | |
void | tlapack::larf (side_t side, direction_t direction, storage_t storeMode, vector_t const &v, const tau_t &tau, matrix_t &C) |
Applies an elementary reflector H to a m-by-n matrix C. | |
template<TLAPACK_SIDE side_t, TLAPACK_STOREV storage_t, TLAPACK_VECTOR vector_t, TLAPACK_SCALAR tau_t, TLAPACK_VECTOR vectorC0_t, TLAPACK_MATRIX matrixC1_t, enable_if_t< std::is_convertible_v< storage_t, StoreV >, int > = 0> | |
void | tlapack::larf (side_t side, storage_t storeMode, vector_t const &x, const tau_t &tau, vectorC0_t &C0, matrixC1_t &C1) |
Applies an elementary reflector defined by tau and v to a m-by-n matrix C decomposed into C0 and C1. | |
template<TLAPACK_SIDE side_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_SVECTOR vector_t, TLAPACK_WORKSPACE work_t, TLAPACK_SCALAR tau_t, TLAPACK_SMATRIX matrix_t, enable_if_t< std::is_convertible_v< direction_t, Direction >, int > = 0> | |
void | tlapack::larf_work (side_t side, direction_t direction, storage_t storeMode, vector_t const &v, const tau_t &tau, matrix_t &C, work_t &work) |
Applies an elementary reflector H to a m-by-n matrix C. Workspace is provided as an argument. | |
template<TLAPACK_SIDE side_t, TLAPACK_STOREV storage_t, TLAPACK_VECTOR vector_t, TLAPACK_WORKSPACE work_t, TLAPACK_SCALAR tau_t, TLAPACK_VECTOR vectorC0_t, TLAPACK_MATRIX matrixC1_t, enable_if_t< std::is_convertible_v< storage_t, StoreV >, int > = 0> | |
void | tlapack::larf_work (side_t side, storage_t storeMode, vector_t const &x, const tau_t &tau, vectorC0_t &C0, matrixC1_t &C1, work_t &work) |
Applies an elementary reflector defined by tau and v to a m-by-n matrix C decomposed into C0 and C1. Workspace is provided as an argument. | |
template<TLAPACK_SMATRIX matrixV_t, TLAPACK_MATRIX matrixT_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t> | |
int | tlapack::larfb (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const matrixT_t &Tmatrix, matrixC_t &C) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right. | |
template<TLAPACK_SMATRIX matrixV_t, TLAPACK_MATRIX matrixT_t, TLAPACK_SMATRIX matrixC_t, TLAPACK_WORKSPACE work_t, TLAPACK_SIDE side_t, TLAPACK_OP trans_t, TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t> | |
int | tlapack::larfb_work (side_t side, trans_t trans, direction_t direction, storage_t storeMode, const matrixV_t &V, const matrixT_t &Tmatrix, matrixC_t &C, work_t &work) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right. Workspace is provided as an argument. | |
template<TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_VECTOR vector_t, enable_if_t< std::is_convertible_v< direction_t, Direction >, int > = 0> | |
void | tlapack::larfg (direction_t direction, storage_t storeMode, vector_t &v, type_t< vector_t > &tau) |
Generates a elementary Householder reflection. | |
template<TLAPACK_STOREV storage_t, TLAPACK_VECTOR vector_t> | |
void | tlapack::larfg (storage_t storeMode, type_t< vector_t > &alpha, vector_t &x, type_t< vector_t > &tau) |
Generates a elementary Householder reflection. | |
template<TLAPACK_DIRECTION direction_t, TLAPACK_STOREV storage_t, TLAPACK_SMATRIX matrixV_t, TLAPACK_VECTOR vector_t, TLAPACK_SMATRIX matrixT_t> | |
int | tlapack::larft (direction_t direction, storage_t storeMode, const matrixV_t &V, const vector_t &tau, matrixT_t &T) |
Forms the triangular factor T of a block reflector H of order n, which is defined as a product of k elementary reflectors. | |
template<int idist, TLAPACK_VECTOR vector_t, class iseed_t > | |
void | tlapack::larnv (iseed_t &iseed, vector_t &x) |
Returns a vector of n random numbers from a uniform or normal distribution. | |
template<TLAPACK_MATRIX matrix_t, TLAPACK_REAL a_type, TLAPACK_REAL b_type, enable_if_t<(is_real< a_type > &&is_real< b_type >), int > = 0> | |
int | tlapack::lascl (BandAccess accessType, const b_type &b, const a_type &a, matrix_t &A) |
Multiplies a matrix A by the real scalar a/b. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t, TLAPACK_REAL a_type, TLAPACK_REAL b_type, enable_if_t<(is_real< a_type > &&is_real< b_type >), int > = 0> | |
int | tlapack::lascl (uplo_t uplo, const b_type &b, const a_type &a, matrix_t &A) |
Multiplies a matrix A by the real scalar a/b. | |
template<TLAPACK_UPLO uplo_t, TLAPACK_MATRIX matrix_t> | |
void | tlapack::laset (uplo_t uplo, const type_t< matrix_t > &alpha, const type_t< matrix_t > &beta, matrix_t &A) |
Initializes a matrix to diagonal and off-diagonal values. | |
template<TLAPACK_VECTOR vector_t> | |
void | tlapack::lassq (const vector_t &x, real_type< type_t< vector_t > > &scale, real_type< type_t< vector_t > > &sumsq) |
Updates a sum of squares represented in scaled form. | |
template<class abs_f , TLAPACK_VECTOR vector_t> | |
void | tlapack::lassq (const vector_t &x, real_type< type_t< vector_t > > &scale, real_type< type_t< vector_t > > &sumsq, abs_f absF) |
Updates a sum of squares represented in scaled form. | |
template<TLAPACK_MATRIX matrixX_t, TLAPACK_MATRIX matrixT_t, TLAPACK_MATRIX matrixB_t, enable_if_t< is_real< type_t< matrixX_t > > &&is_real< type_t< matrixT_t > > &&is_real< type_t< matrixB_t > >, bool > = true> | |
int | tlapack::lasy2 (Op trans_l, Op trans_r, int isign, const matrixT_t &TL, const matrixT_t &TR, const matrixB_t &B, type_t< matrixX_t > &scale, matrixX_t &X, type_t< matrixX_t > &xnorm) |
lasy2 solves the Sylvester matrix equation where the matrices are of order 1 or 2. | |
template<TLAPACK_SMATRIX matrix_t> | |
void | tlapack::lu_mult (matrix_t &A, const LuMultOpts &opts={}) |
in-place multiplication of lower triangular matrix L and upper triangular matrix U. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_CVECTOR vector_t> | |
void | tlapack::move_bulge (matrix_t &H, vector_t &v, complex_type< type_t< matrix_t > > s1, complex_type< type_t< matrix_t > > s2) |
Given a 4-by-3 matrix H and small order reflector v, move_bulge applies the delayed right update to the last row and calculates a new reflector to move the bulge down. | |
template<TLAPACK_SMATRIX matrix_t, TLAPACK_SVECTOR alpha_t, TLAPACK_SVECTOR beta_t> | |
int | tlapack::multishift_qz (bool want_s, bool want_q, bool want_z, size_type< matrix_t > ilo, size_type< matrix_t > ihi, matrix_t &A, matrix_t &B, alpha_t &alpha, beta_t &beta, matrix_t &Q, matrix_t &Z, FrancisOpts &opts) |
multishift_qz computes the eigenvalues of a matrix pair (H,T), where H is an upper Hessenberg matrix and T is upper triangular, using the multishift QZ algorithm with AED. | |
template<TLAPACK_VECTOR vector_t, TLAPACK_REAL alpha_t, enable_if_t< is_real< alpha_t >, int > = 0> | |
void | tlapack::rscl (const alpha_t &alpha, vector_t &x) |
Scale vector by the reciprocal of a constant, \(x := x / \alpha\). | |
template<TLAPACK_MATRIX matrix_t> | |
int | tlapack::schur_move (bool want_q, matrix_t &A, matrix_t &Q, size_type< matrix_t > &ifst, size_type< matrix_t > &ilst) |
schur_move reorders the Schur factorization of a matrix S = Q*A*Q**H, so that the diagonal element of S with row index IFST is moved to row ILST. | |
template<TLAPACK_CSMATRIX matrix_t, enable_if_t< is_real< type_t< matrix_t > >, bool > = true> | |
int | tlapack::schur_swap (bool want_q, matrix_t &A, matrix_t &Q, const size_type< matrix_t > &j0, const size_type< matrix_t > &n1, const size_type< matrix_t > &n2) |
schur_swap, swaps 2 eigenvalues of A. | |
template<typename T , enable_if_t<!is_complex< T >, bool > = true> | |
void | tlapack::singularvalues22 (const T &f, const T &g, const T &h, T &ssmin, T &ssmax) |
Computes the singular value decomposition of a 2-by-2 real triangular matrix. | |
template<typename T , enable_if_t<!is_complex< T >, bool > = true> | |
void | tlapack::svd22 (const T &f, const T &g, const T &h, T &ssmin, T &ssmax, T &csl, T &snl, T &csr, T &snr) |
Computes the singular value decomposition of a 2-by-2 real triangular matrix. | |
template<TLAPACK_SMATRIX matrixA_t, TLAPACK_SMATRIX matrixB_t> | |
void | tlapack::transpose (matrixA_t &A, matrixB_t &B, const TransposeOpts &opts={}) |
transpose a matrix A into a matrix B. | |
template<TLAPACK_SMATRIX matrix_t> | |
int | tlapack::ul_mult (matrix_t &A) |
ul_mult computes the matrix product of an upper triangular matrix U and a lower triangular unital matrix L Given input matrix A, nonzero part of L is the subdiagonal of A and on the diagonal of L is assumed to be 1, and the nonzero part of U is diagonal and super-diagonal part of A | |
Auxiliary routines that are used by the computational routines.
real_type< type_t< matrix_t > > tlapack::check_generalized_similarity_transform | ( | matrix_t & | A, |
matrix_t & | Q, | ||
matrix_t & | Z, | ||
matrix_t & | B | ||
) |
Calculates ||Q'*A*Z - B||.
[in] | A | n by n matrix |
[in] | Q | n by n unitary matrix |
[in] | Z | n by n unitary matrix |
[in] | B | n by n matrix |
real_type< type_t< matrix_t > > tlapack::check_generalized_similarity_transform | ( | matrix_t & | A, |
matrix_t & | Q, | ||
matrix_t & | Z, | ||
matrix_t & | B, | ||
matrix_t & | res, | ||
matrix_t & | work | ||
) |
Calculates res = Q'*A*Z - B and the frobenius norm of res.
[in] | A | n by n matrix |
[in] | Q | n by n unitary matrix |
[in] | Z | n by n unitary matrix |
[in] | B | n by n matrix |
[out] | res | n by n matrix as defined above |
[out] | work | n by n workspace matrix |
real_type< type_t< matrix_t > > tlapack::check_orthogonality | ( | matrix_t & | Q | ) |
Calculates ||Q'*Q - I||_F if m <= n or ||Q*Q' - I||_F otherwise.
[in] | Q | m by n (almost) orthogonal matrix |
real_type< type_t< matrix_t > > tlapack::check_orthogonality | ( | matrix_t & | Q, |
matrix_t & | res | ||
) |
Calculates res = Q'*Q - I if m <= n or res = Q*Q' otherwise Also computes the frobenius norm of res.
[in] | Q | m by n (almost) orthogonal matrix |
[out] | res | n by n matrix as defined above |
real_type< type_t< matrix_t > > tlapack::check_similarity_transform | ( | matrix_t & | A, |
matrix_t & | Q, | ||
matrix_t & | B | ||
) |
Calculates ||Q'*A*Q - B||.
[in] | A | n by n matrix |
[in] | Q | n by n unitary matrix |
[in] | B | n by n matrix |
real_type< type_t< matrix_t > > tlapack::check_similarity_transform | ( | matrix_t & | A, |
matrix_t & | Q, | ||
matrix_t & | B, | ||
matrix_t & | res, | ||
matrix_t & | work | ||
) |
Calculates res = Q'*A*Q - B and the frobenius norm of res.
[in] | A | n by n matrix |
[in] | Q | n by n unitary matrix |
[in] | B | n by n matrix |
[out] | res | n by n matrix as defined above |
[out] | work | n by n workspace matrix |
void tlapack::conjtranspose | ( | matrixA_t & | A, |
matrixB_t & | B, | ||
const TransposeOpts & | opts = {} |
||
) |
conjugate transpose a matrix A into a matrix B.
[in] | A | m-by-n matrix The matrix to be transposed |
[out] | B | n-by-m matrix On exit, B = A**H |
[in] | opts | Options. |
void tlapack::conjugate | ( | vector_t & | x | ) |
Conjugates a vector.
[in] | x | Vector of size n. |
int tlapack::generalized_schur_move | ( | bool | want_q, |
bool | want_z, | ||
matrix_t & | A, | ||
matrix_t & | B, | ||
matrix_t & | Q, | ||
matrix_t & | Z, | ||
size_type< matrix_t > & | ifst, | ||
size_type< matrix_t > & | ilst | ||
) |
generalized_schur_move reorders the generalized Schur factorization of a pencil ( S, T ) = Q(A,B)Z**H so that the diagonal elements of (S,T) with row index IFST are moved to row ILST.
[in] | want_q | bool Whether or not to apply the transformations to Q |
[in] | want_z | bool Whether or not to apply the transformations to Q |
[in,out] | A | n-by-n matrix. Must be in Schur form |
[in,out] | B | n-by-n matrix. Must be in Schur form |
[in,out] | Q | n-by-n matrix. unitary matrix, not referenced if want_q is false |
[in,out] | Z | n-by-n matrix. unitary matrix, not referenced if want_q is false |
[in,out] | ifst | integer Initial row index of the eigenvalue block |
[in,out] | ilst | integer Index of the eigenvalue block at the exit of the routine. If it is not possible to move the eigenvalue block to the given ilst, the value may be changed. |
int tlapack::generalized_schur_swap | ( | bool | want_q, |
bool | want_z, | ||
matrix_t & | A, | ||
matrix_t & | B, | ||
matrix_t & | Q, | ||
matrix_t & | Z, | ||
const size_type< matrix_t > & | j0, | ||
const size_type< matrix_t > & | n1, | ||
const size_type< matrix_t > & | n2 | ||
) |
schur_swap, swaps 2 eigenvalues of A.
[in] | want_q | bool Whether or not to apply the transformations to Q |
[in] | want_z | bool Whether or not to apply the transformations to Z |
[in,out] | A | n-by-n matrix. Must be in Schur form |
[in,out] | B | n-by-n matrix. Must be in Schur form |
[in,out] | Q | n-by-n matrix. Orthogonal matrix, not referenced if want_q is false |
[in,out] | Z | n-by-n matrix. Orthogonal matrix, not referenced if want_q is false |
[in] | j0 | integer Index of first eigenvalue block |
[in] | n1 | integer Size of first eigenvalue block |
[in] | n2 | integer Size of second eigenvalue block |
Implementation for complex matrices
auto tlapack::infnorm_colmajor | ( | const matrix_t & | A | ) |
Calculates the infinity norm of a column-major matrix.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | A | m-by-n matrix. |
auto tlapack::infnorm_colmajor_work | ( | const matrix_t & | A, |
work_t & | work | ||
) |
Calculates the infinity norm of a column-major matrix. Workspace is provided as an argument.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | A | m-by-n matrix. |
work | Workspace. Use the workspace query to determine the size needed. |
auto tlapack::infnorm_hermitian_colmajor | ( | uplo_t | uplo, |
const matrix_t & | A | ||
) |
Calculates the infinity norm of a column-major hermitian matrix.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | A | n-by-n matrix. |
Calculates the infinity norm of a column-major hermitian matrix. Workspace is provided as an argument.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | A | n-by-n matrix. |
work | Workspace. Use the workspace query to determine the size needed. |
auto tlapack::infnorm_symmetric_colmajor | ( | uplo_t | uplo, |
const matrix_t & | A | ||
) |
Calculates the infinity norm of a column-major symmetric matrix.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | A | n-by-n matrix. |
Calculates the infinity norm of a column-major symmetric matrix. Workspace is provided as an argument.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | A | n-by-n matrix. |
work | Workspace. Use the workspace query to determine the size needed. |
Calculates the infinity norm of a column-major triangular matrix.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | A | m-by-n matrix. |
auto tlapack::infnorm_triangular_colmajor_work | ( | uplo_t | uplo, |
diag_t | diag, | ||
const matrix_t & | A, | ||
work_t & | work | ||
) |
Calculates the infinity norm of a column-major triangular matrix. Workspace is provided as an argument.
Code optimized for the infinity norm on column-major layouts using a workspace of size at least m, where m is the number of rows of A.
[in] | uplo |
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | A | m-by-n matrix. |
work | Workspace. Use the workspace query to determine the size needed. |
void tlapack::inv_house3 | ( | const matrix_t & | A, |
vector_t & | v, | ||
type_t< vector_t > & | tau | ||
) |
Inv_house calculates a reflector to reduce the first column in a 2x3 matrix A from the right to zero.
Note that this is special because a reflector applied from the right would usually reduce a row, not a column. This is known as an inverse reflector.
We need to solve the system A x = 0 If we then calculate a reflector that reduces x: H x = alpha e_1, then H is the reflector that we were looking for.
Because the reflector is invariant w.r.t. the scale of x, we will solve a scaled system so that x0 = 1. The rest of x can then be calculated using: [A11 A12] [x1] = - scale * [A10] [A21 A22] [x2] [A20]
We solve this system of equations using fully pivoted LU
[in] | A | 2x3 matrix. |
[out] | v | vector of size 3 |
[out] | tau |
int tlapack::labrd | ( | A_t & | A, |
vector_t & | tauq, | ||
vector_t & | taup, | ||
X_t & | X, | ||
matrixY_t & | Y | ||
) |
Reduces the first nb rows and columns of a general m by n matrix A to upper or lower bidiagonal form by an unitary transformation Q**H * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A.
If m >= n, A is reduced to upper bidiagonal form; if m < n, to lower bidiagonal form.
This is an auxiliary routine called by gebrd
[in,out] | A | m-by-n matrix. |
[out] | tauq | vector of length nb. The scalar factors of the elementary reflectors which represent the unitary matrix Q. |
[out] | taup | vector of length nb. The scalar factors of the elementary reflectors which represent the unitary matrix P. |
[out] | X | m-by-nb matrix. |
[out] | Y | n-by-nb matrix. |
Copies a matrix from A to B.
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
A | m-by-n matrix. | |
B | matrix with at least m rows and at least n columns. |
void tlapack::ladiv | ( | const real_t & | a, |
const real_t & | b, | ||
const real_t & | c, | ||
const real_t & | d, | ||
real_t & | p, | ||
real_t & | q | ||
) |
Performs complex division in real arithmetic.
\[ p + iq = (a + ib) / (c + id) \]
real_t | Floating-point type. |
[in] | a | Real part of numerator. |
[in] | b | Imaginary part of numerator. |
[in] | c | Real part of denominator. |
[in] | d | Imaginary part of denominator. |
[out] | p | Real part of quotient. |
[out] | q | Imaginary part of quotient. |
Performs complex division in real arithmetic with complex arguments.
real_t | Floating-point type. |
[in] | x | Complex numerator. |
[in] | y | Complex denominator. |
int tlapack::lahqr | ( | 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 | ||
) |
lahqr computes the eigenvalues and optionally the Schur factorization of an upper Hessenberg matrix, using the double-shift implicit QR algorithm.
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. |
Implementation for complex matrices.
void tlapack::lahqr_eig22 | ( | T | a00, |
T | a01, | ||
T | a10, | ||
T | a11, | ||
complex_type< T > & | s1, | ||
complex_type< T > & | s2 | ||
) |
Computes the eigenvalues of a 2x2 matrix A.
[in] | a00 | Element (0,0) of A. |
[in] | a01 | Element (0,1) of A. |
[in] | a10 | Element (1,0) of A. |
[in] | a11 | Element (1,1) of A. |
[out] | s1 | |
[out] | s2 | s1 and s2 are the eigenvalues of A |
void tlapack::lahqr_schur22 | ( | T & | a, |
T & | b, | ||
T & | c, | ||
T & | d, | ||
complex_type< T > & | s1, | ||
complex_type< T > & | s2, | ||
T & | cs, | ||
T & | sn | ||
) |
Computes the Schur factorization of a 2x2 matrix A.
A = [a b] = [cs -sn] [aa bb] [ cs sn] [c d] [sn cs] [cc dd] [-sn cs]
[in,out] | a | scalar, A(0,0). |
[in,out] | b | scalar, A(0,1). |
[in,out] | c | scalar, A(1,0). |
[in,out] | d | scalar, A(1,1). On entry, the elements of the matrix A. On exit, the elements of the Schur factor (aa, bb, cc and dd). |
[out] | s1 | complex scalar. First eigenvalue of the matrix. |
[out] | s2 | complex scalar. Second eigenvalue of the matrix. |
[out] | cs | scalar. Cosine factor of the rotation |
[out] | sn | scalar. Sine factor of the rotation |
int tlapack::lahqr_shiftcolumn | ( | const matrix_t & | H, |
vector_t & | v, | ||
complex_type< type_t< matrix_t > > | s1, | ||
complex_type< type_t< matrix_t > > | s2 | ||
) |
Given a 2-by-2 or 3-by-3 matrix H, lahqr_shiftcolumn calculates a multiple of the product: (H - s1*I)*(H - s2*I)*e1.
This is used to introduce shifts in the QR algorithm
[in] | H | 2x2 or 3x3 matrix. The matrix H as in the formula above. |
[out] | v | vector of size 2 or 3 On exit, a multiple of the product |
[in] | s1 | The scalar s1 as in the formula above |
[in] | s2 | The scalar s2 as in the formula above |
int tlapack::lahqr_shiftcolumn | ( | const matrix_t & | H, |
vector_t & | v, | ||
type_t< matrix_t > | s1, | ||
type_t< matrix_t > | s2 | ||
) |
Given a 2-by-2 or 3-by-3 matrix H, lahqr_shiftcolumn calculates a multiple of the product: (H - s1*I)*(H - s2*I)*e1.
This is used to introduce shifts in the QR algorithm
[in] | H | 2x2 or 3x3 matrix. The matrix H as in the formula above. |
[out] | v | vector of size 2 or 3 On exit, a multiple of the product |
[in] | s1 | The scalar s1 as in the formula above |
[in] | s2 | The scalar s2 as in the formula above |
int tlapack::lahqz | ( | bool | want_s, |
bool | want_q, | ||
bool | want_z, | ||
size_type< matrix_t > | ilo, | ||
size_type< matrix_t > | ihi, | ||
matrix_t & | A, | ||
matrix_t & | B, | ||
alpha_t & | alpha, | ||
beta_t & | beta, | ||
matrix_t & | Q, | ||
matrix_t & | Z | ||
) |
lahqz computes the eigenvalues of a matrix pair (H,T), where H is an upper Hessenberg matrix and T is upper triangular, using the single/double-shift implicit QZ algorithm.
[in] | want_s | bool. If true, the full Schur form 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. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n. |
[in,out] | A | n by n matrix. |
[in,out] | B | n by n matrix. |
[out] | alpha | size n vector. |
[out] | beta | size n vector. |
[in,out] | Q | n by n matrix. |
[in,out] | Z | n by n matrix. |
int tlapack::lahqz_shiftcolumn | ( | const matrix_t & | A, |
const matrix_t & | B, | ||
vector_t & | v, | ||
complex_type< type_t< matrix_t > > | s1, | ||
complex_type< type_t< matrix_t > > | s2, | ||
type_t< matrix_t > | beta1, | ||
type_t< matrix_t > | beta2 | ||
) |
Given a 2-by-2 or 3-by-3 matrix pencil (A,B), lahqz_shiftcolumn calculates a multiple of the product: (beta2*A - s2*B)*B^(-1)*(beta1*A - s1*B)*B^(-1)*e1.
This is used to introduce shifts in the QZ algorithm
[in] | A | 2x2 or 3x3 matrix. The matrix A as in the formula above. |
[in] | B | 2x2 or 3x3 matrix. The matrix B as in the formula above. |
[out] | v | vector of size 2 or 3 On exit, a multiple of the product |
[in] | s1 | The scalar s1 as in the formula above |
[in] | s2 | The scalar s2 as in the formula above |
[in] | beta1 | The scalar beta1 as in the formula above |
[in] | beta2 | The scalar beta2 as in the formula above |
int tlapack::lahr2 | ( | size_type< matrix_t > | k, |
size_type< matrix_t > | nb, | ||
matrix_t & | A, | ||
vector_t & | tau, | ||
matrixT_t & | T, | ||
matrixY_t & | Y | ||
) |
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].
[in] | k | integer |
[in] | nb | integer |
[in,out] | A | n-by-n matrix. |
[out] | tau | Real vector of length n-1. The scalar factors of the elementary reflectors. |
[out] | T | nb-by-nb matrix. |
[out] | Y | n-by-nb matrix. |
auto tlapack::lange | ( | norm_t | normType, |
const matrix_t & | A | ||
) |
Calculates the norm of a matrix.
norm_t | Either Norm or any class that implements operator Norm() . |
[in] | normType |
|
[in] | A | m-by-n matrix. |
Calculates the norm of a hermitian matrix.
norm_t | Either Norm or any class that implements operator Norm() . |
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | normType |
|
[in] | uplo |
|
[in] | A | n-by-n hermitian matrix. |
Calculates the norm of a symmetric matrix.
norm_t | Either Norm or any class that implements operator Norm() . |
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | normType |
|
[in] | uplo |
|
[in] | A | n-by-n symmetric matrix. |
Calculates the norm of a symmetric matrix.
norm_t | Either Norm or any class that implements operator Norm() . |
uplo_t | Either Uplo or any class that implements operator Uplo() . |
diag_t | Either Diag or any class that implements operator Diag() . |
[in] | normType |
|
[in] | uplo |
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | A | m-by-n triangular matrix. |
Finds \(\sqrt{x^2+y^2}\), taking care not to cause unnecessary overflow.
[in] | x | scalar value x |
[in] | y | scalar value y |
real_type< TX, TY, TZ > tlapack::lapy3 | ( | const TX & | x, |
const TY & | y, | ||
const TZ & | z | ||
) |
Finds \(\sqrt{x^2+y^2+z^2}\), taking care not to cause unnecessary overflow or unnecessary underflow.
real_t | Floating-point type. |
[in] | x | scalar value x |
[in] | y | scalar value y |
[in] | z | scalar value y |
void tlapack::larf | ( | side_t | side, |
direction_t | direction, | ||
storage_t | storeMode, | ||
vector_t const & | v, | ||
const tau_t & | tau, | ||
matrix_t & | C | ||
) |
Applies an elementary reflector H to a m-by-n matrix C.
The elementary reflector H can be applied on either the left or right, with
\[ H = I - \tau v v^H. \]
where v = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward.
side_t | Either Side or any class that implements operator Side() . |
[in] | side |
|
[in] | direction | v = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward. |
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in] | v | Vector of size m if side = Side::Left, or n if side = Side::Right. |
[in] | tau | Value of tau in the representation of H. |
[in,out] | C | On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) if side = Side::Left, or \(C H\) if side = Side::Right. |
void tlapack::larf | ( | side_t | side, |
storage_t | storeMode, | ||
vector_t const & | x, | ||
const tau_t & | tau, | ||
vectorC0_t & | C0, | ||
matrixC1_t & | C1 | ||
) |
Applies an elementary reflector defined by tau and v to a m-by-n matrix C decomposed into C0 and C1.
\[ C_0 = (1-\tau) C_0 - \tau v^H C_1, \\ C_1 = -\tau v C_0 + (I-\tau vv^H) C_1, \]
if side = Side::Left, or
\[ C_0 = (1-\tau) C_0 -\tau C_1 v, \\ C_1 = -\tau C_0 v^H + C_1 (I-\tau vv^H), \]
if side = Side::Right.
The elementary reflector is defined as
\[ H = \begin{bmatrix} 1-\tau & -\tau v^H \\ -\tau v & I-\tau vv^H \end{bmatrix} \]
side_t | Either Side or any class that implements operator Side() . |
[in] | side |
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in] | x | Vector \(v\) if storeMode = StoreV::Columnwise, or \(v^H\) if storeMode = StoreV::Rowwise. |
[in] | tau | Value of tau in the representation of H. |
[in,out] | C0 | Vector of size n if side = Side::Left, or m if side = Side::Right. On exit, C0 is overwritten by
|
[in,out] | C1 | Matrix of size (m-1)-by-n if side = Side::Left, or m-by-(n-1) if side = Side::Right. On exit, C1 is overwritten by
|
void tlapack::larf_work | ( | side_t | side, |
direction_t | direction, | ||
storage_t | storeMode, | ||
vector_t const & | v, | ||
const tau_t & | tau, | ||
matrix_t & | C, | ||
work_t & | work | ||
) |
Applies an elementary reflector H to a m-by-n matrix C. Workspace is provided as an argument.
The elementary reflector H can be applied on either the left or right, with
\[ H = I - \tau v v^H. \]
where v = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward.
side_t | Either Side or any class that implements operator Side() . |
[in] | side |
|
[in] | direction | v = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward. |
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in] | v | Vector of size m if side = Side::Left, or n if side = Side::Right. |
[in] | tau | Value of tau in the representation of H. |
[in,out] | C | On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) if side = Side::Left, or \(C H\) if side = Side::Right. |
work | Workspace. Use the workspace query to determine the size needed. |
void tlapack::larf_work | ( | side_t | side, |
storage_t | storeMode, | ||
vector_t const & | x, | ||
const tau_t & | tau, | ||
vectorC0_t & | C0, | ||
matrixC1_t & | C1, | ||
work_t & | work | ||
) |
Applies an elementary reflector defined by tau and v to a m-by-n matrix C decomposed into C0 and C1. Workspace is provided as an argument.
\[ C_0 = (1-\tau) C_0 - \tau v^H C_1, \\ C_1 = -\tau v C_0 + (I-\tau vv^H) C_1, \]
if side = Side::Left, or
\[ C_0 = (1-\tau) C_0 -\tau C_1 v, \\ C_1 = -\tau C_0 v^H + C_1 (I-\tau vv^H), \]
if side = Side::Right.
The elementary reflector is defined as
\[ H = \begin{bmatrix} 1-\tau & -\tau v^H \\ -\tau v & I-\tau vv^H \end{bmatrix} \]
side_t | Either Side or any class that implements operator Side() . |
[in] | side |
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in] | x | Vector \(v\) if storeMode = StoreV::Columnwise, or \(v^H\) if storeMode = StoreV::Rowwise. |
[in] | tau | Value of tau in the representation of H. |
[in,out] | C0 | Vector of size n if side = Side::Left, or m if side = Side::Right. On exit, C0 is overwritten by
|
[in,out] | C1 | Matrix of size (m-1)-by-n if side = Side::Left, or m-by-(n-1) if side = Side::Right. On exit, C1 is overwritten by
|
work | Workspace. Use the workspace query to determine the size needed. |
int tlapack::larfb | ( | side_t | side, |
trans_t | trans, | ||
direction_t | direction, | ||
storage_t | storeMode, | ||
const matrixV_t & | V, | ||
const matrixT_t & | Tmatrix, | ||
matrixC_t & | C | ||
) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right.
[in] | side |
|
[in] | trans |
|
[in] | direction | Indicates how H is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored:
|
[in] | V |
|
[in] | Tmatrix | The k-by-k matrix T. The triangular k-by-k matrix T in the representation of the block reflector. |
[in,out] | C | On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) or \(H^H C\) or \(C H\) or \(C H^H\). |
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 n = 5 and k = 3. The elements equal to 1 are not stored. The rest of the array 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 )
int tlapack::larfb_work | ( | side_t | side, |
trans_t | trans, | ||
direction_t | direction, | ||
storage_t | storeMode, | ||
const matrixV_t & | V, | ||
const matrixT_t & | Tmatrix, | ||
matrixC_t & | C, | ||
work_t & | work | ||
) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right. Workspace is provided as an argument.
[in] | side |
|
[in] | trans |
|
[in] | direction | Indicates how H is formed from a product of elementary reflectors.
|
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored:
|
[in] | V |
|
[in] | Tmatrix | The k-by-k matrix T. The triangular k-by-k matrix T in the representation of the block reflector. |
[in,out] | C | On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) or \(H^H C\) or \(C H\) or \(C H^H\). |
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 n = 5 and k = 3. The elements equal to 1 are not stored. The rest of the array 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. |
void tlapack::larfg | ( | direction_t | direction, |
storage_t | storeMode, | ||
vector_t & | v, | ||
type_t< vector_t > & | tau | ||
) |
Generates a elementary Householder reflection.
larfg generates a elementary Householder reflection H of order n, such that
H' * ( alpha ) = ( beta ), H' * H = I. ( x ) ( 0 )
if storeMode == StoreV::Columnwise, or
( alpha x ) * H = ( beta 0 ), H' * H = I.
if storeMode == StoreV::Rowwise, where alpha and beta are scalars, with beta real, and x is an (n-1)-element vector. H is represented in the form
H = I - tau * ( 1 ) * ( 1 y' ) ( y )
if storeMode == StoreV::Columnwise, or
H = I - tau * ( 1 ) * ( 1 y ) ( y' )
if storeMode == StoreV::Rowwise, where tau is a scalar and y is a (n-1)-element vector. Note that H is symmetric but not hermitian.
If the elements of x are all zero and alpha is real, then tau = 0 and H is taken to be the identity matrix.
Otherwise 1 <= real(tau) <= 2 and abs(tau-1) <= 1.
direction_t | Either Direction or any class that implements operator Direction() . |
storage_t | Either StoreV or any class that implements operator StoreV() . |
[in] | direction | v = [ alpha x ] if direction == Direction::Forward and v = [ x alpha ] if direction == Direction::Backward. |
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in,out] | v | Vector of length n. On entry, v = [ alpha x ] if direction == Direction::Forward and v = [ x alpha ] if direction == Direction::Backward. On exit, v = [ 1 y ] if direction == Direction::Forward and v = [ y 1 ] if direction == Direction::Backward. |
[out] | tau | The value tau. |
void tlapack::larfg | ( | storage_t | storeMode, |
type_t< vector_t > & | alpha, | ||
vector_t & | x, | ||
type_t< vector_t > & | tau | ||
) |
Generates a elementary Householder reflection.
larfg generates a elementary Householder reflection H of order n, such that
H' * ( alpha ) = ( beta ), H' * H = I. ( x ) ( 0 )
if storeMode == StoreV::Columnwise, or
( alpha x ) * H = ( beta 0 ), H' * H = I.
if storeMode == StoreV::Rowwise, where alpha and beta are scalars, with beta real, and x is an (n-1)-element vector. H is represented in the form
H = I - tau * ( 1 ) * ( 1 y' ) ( y )
if storeMode == StoreV::Columnwise, or
H = I - tau * ( 1 ) * ( 1 y ) ( y' )
if storeMode == StoreV::Rowwise, where tau is a scalar and y is a (n-1)-element vector. Note that H is symmetric but not hermitian.
If the elements of x are all zero and alpha is real, then tau = 0 and H is taken to be the identity matrix.
Otherwise 1 <= real(tau) <= 2 and abs(tau-1) <= 1.
[in] | storeMode | Indicates how the vectors which define the elementary reflectors are stored: |
[in,out] | alpha | On entry, the value alpha. On exit, it is overwritten with the value beta. |
[in,out] | x | Vector of length n-1. On entry, the vector x. On exit, it is overwritten with the vector y. |
[out] | tau | The value tau. |
int tlapack::larft | ( | direction_t | direction, |
storage_t | storeMode, | ||
const matrixV_t & | V, | ||
const vector_t & | tau, | ||
matrixT_t & | T | ||
) |
Forms the triangular factor T of a block reflector H of order n, which is defined as a product of k elementary reflectors.
If direction = Direction::Forward, H = H_1 H_2 ... H_k and T is upper triangular. If direction = Direction::Backward, H = H_k ... H_2 H_1 and T is lower triangular.
If storeMode = StoreV::Columnwise, the vector which defines the elementary reflector H(i) is stored in the i-th column of the array V, and
H = I - V * T * V'
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 n = 5 and k = 3. The elements equal to 1 are not stored. The rest of the array 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 )
direction_t | Either Direction or any class that implements operator Direction() . |
storage_t | Either StoreV or any class that implements operator StoreV() . |
[in] | direction | Indicates how H 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 containing the scalar factors of the elementary reflectors H. |
[out] | T | Matrix of size k-by-k containing the triangular factors of the block reflector.
|
Returns a vector of n random numbers from a uniform or normal distribution.
This implementation uses the Mersenne Twister 19937 generator (std::mt19937), which is a Mersenne Twister pseudo-random generator of 32-bit numbers with a state size of 19937 bits.
Requires ISO C++ 2011 random number generators.
idist | Specifies the distribution:
|
[in,out] | iseed | Seed for the random number generator. The seed is updated inside the routine ( seed := seed + 1 ). |
[out] | x | Vector of length n. |
int tlapack::lascl | ( | BandAccess | accessType, |
const b_type & | b, | ||
const a_type & | a, | ||
matrix_t & | A | ||
) |
Multiplies a matrix A by the real scalar a/b.
Specific implementation for band access types.
[in] | accessType | Determines the entries of A that are scaled by a/b. |
[in] | b | The denominator of the scalar a/b. |
[in] | a | The numerator of the scalar a/b. |
[in,out] | A | Matrix to be scaled by a/b. |
int tlapack::lascl | ( | uplo_t | uplo, |
const b_type & | b, | ||
const a_type & | a, | ||
matrix_t & | A | ||
) |
Multiplies a matrix A by the real scalar a/b.
Multiplication of a matrix A by scalar a/b is done without over/underflow as long as the final result \(a A/b\) does not over/underflow.
uplo_t | Type of access inside the algorithm. Either Uplo or any type that implements operator Uplo(). |
matrix_t | Matrix type. |
a_type | Type of the coefficient a. a_type cannot be a complex type. |
b_type | Type of the coefficient b. b_type cannot be a complex type. |
[in] | uplo | Determines the entries of A that are scaled by a/b. The following access types are allowed: Uplo::General, Uplo::UpperHessenberg, Uplo::LowerHessenberg, Uplo::Upper, Uplo::Lower, Uplo::StrictUpper, Uplo::StrictLower. |
[in] | b | The denominator of the scalar a/b. |
[in] | a | The numerator of the scalar a/b. |
[in,out] | A | Matrix to be scaled by a/b. |
void tlapack::laset | ( | uplo_t | uplo, |
const type_t< matrix_t > & | alpha, | ||
const type_t< matrix_t > & | beta, | ||
matrix_t & | A | ||
) |
Initializes a matrix to diagonal and off-diagonal values.
uplo_t | Either Uplo or any class that implements operator Uplo() . |
[in] | uplo |
|
[in] | alpha | Value to assign to the off-diagonal elements of A. |
[in] | beta | Value to assign to the diagonal elements of A. |
[out] | A | m-by-n matrix. |
void tlapack::lassq | ( | const vector_t & | x, |
real_type< type_t< vector_t > > & | scale, | ||
real_type< type_t< vector_t > > & | sumsq | ||
) |
Updates a sum of squares represented in scaled form.
\[ scl smsq := \sum_{i = 0}^n |x_i|^2 + scale^2 sumsq. \]
Specific implementation using
where T is the type_t< vector_t >.
void tlapack::lassq | ( | const vector_t & | x, |
real_type< type_t< vector_t > > & | scale, | ||
real_type< type_t< vector_t > > & | sumsq, | ||
abs_f | absF | ||
) |
Updates a sum of squares represented in scaled form.
\[ scl smsq := \sum_{i = 0}^n |x_i|^2 + scale^2 sumsq, \]
scale and sumsq are assumed to be non-negative.
[in] | x | Vector of size n. |
[in,out] | scale | Real scalar. On entry, the value scale in the equation above. On exit, scale is overwritten with scl , the scaling factor for the sum of squares. |
[in,out] | sumsq | Real scalar. On entry, the value sumsq in the equation above. On exit, sumsq is overwritten with smsq , the basic sum of squares from which scl has been factored out. |
[in] | absF | Lambda function that computes the absolute value. |
int tlapack::lasy2 | ( | Op | trans_l, |
Op | trans_r, | ||
int | isign, | ||
const matrixT_t & | TL, | ||
const matrixT_t & | TR, | ||
const matrixB_t & | B, | ||
type_t< matrixX_t > & | scale, | ||
matrixX_t & | X, | ||
type_t< matrixX_t > & | xnorm | ||
) |
lasy2 solves the Sylvester matrix equation where the matrices are of order 1 or 2.
lasy2 solves for the N1 by N2 matrix X, 1 <= N1,N2 <= 2, in op(TL)*X + ISGN*X*op(TR) = SCALE*B,
where TL is N1 by N1, TR is N2 by N2, B is N1 by N2, and ISGN = 1 or -1. op(T) = T or T**T, where T**T denotes the transpose of T.
scale = 1
and xnorm = 0
if N1 = N2 = 0.void tlapack::lu_mult | ( | matrix_t & | A, |
const LuMultOpts & | opts = {} |
||
) |
in-place multiplication of lower triangular matrix L and upper triangular matrix U.
this is the recursive variant
[in,out] | A | n-by-n matrix On entry, the strictly lower triangular entries of A contain the matrix L. L is assumed to have unit diagonal. The upper triangular entires of A contain the matrix U. On exit, A contains the product L*U. |
[in] | opts | Options. |
void tlapack::move_bulge | ( | matrix_t & | H, |
vector_t & | v, | ||
complex_type< type_t< matrix_t > > | s1, | ||
complex_type< type_t< matrix_t > > | s2 | ||
) |
Given a 4-by-3 matrix H and small order reflector v, move_bulge applies the delayed right update to the last row and calculates a new reflector to move the bulge down.
If the bulge collapses, an attempt is made to reintroduce it using shifts s1 and s2.
[in,out] | H | 4x3 matrix. |
[in,out] | v | vector of size 3 On entry, the delayed reflector to apply The first element of the reflector is assumed to be one, and v[0] instead stores tau. On exit, the reflector that moves the bulge down one position |
[in] | s1 | complex valued shift |
[in] | s2 | complex valued shift |
int tlapack::multishift_qz | ( | bool | want_s, |
bool | want_q, | ||
bool | want_z, | ||
size_type< matrix_t > | ilo, | ||
size_type< matrix_t > | ihi, | ||
matrix_t & | A, | ||
matrix_t & | B, | ||
alpha_t & | alpha, | ||
beta_t & | beta, | ||
matrix_t & | Q, | ||
matrix_t & | Z, | ||
FrancisOpts & | opts | ||
) |
multishift_qz computes the eigenvalues of a matrix pair (H,T), where H is an upper Hessenberg matrix and T is upper triangular, using the multishift QZ algorithm with AED.
[in] | want_s | bool. If true, the full Schur form 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. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n. |
[in,out] | A | n by n matrix. |
[in,out] | B | n by n matrix. |
[out] | alpha | size n vector. |
[out] | beta | size n vector. |
[in,out] | Q | n by n matrix. |
[in,out] | Z | n by n matrix. |
[in,out] | opts | Options. |
Scale vector by the reciprocal of a constant, \(x := x / \alpha\).
If alpha is real, then this routine is equivalent to scal(1/alpha, x). This is done without overflow or underflow as long as the final result x/a does not overflow or underflow.
If alpha is complex, then we use the following algorithm:
[in] | alpha | Scalar. |
[in,out] | x | A n-element vector. |
int tlapack::schur_move | ( | bool | want_q, |
matrix_t & | A, | ||
matrix_t & | Q, | ||
size_type< matrix_t > & | ifst, | ||
size_type< matrix_t > & | ilst | ||
) |
schur_move reorders the Schur factorization of a matrix S = Q*A*Q**H, so that the diagonal element of S with row index IFST is moved to row ILST.
[in] | want_q | bool Whether or not to apply the transformations to Q |
[in,out] | A | n-by-n matrix. Must be in Schur form |
[in,out] | Q | n-by-n matrix. Orthogonal matrix, not referenced if want_q is false |
[in,out] | ifst | integer Initial row index of the eigenvalue block |
[in,out] | ilst | integer Index of the eigenvalue block at the exit of the routine. If it is not possible to move the eigenvalue block to the given ilst, the value may be changed. |
int tlapack::schur_swap | ( | bool | want_q, |
matrix_t & | A, | ||
matrix_t & | Q, | ||
const size_type< matrix_t > & | j0, | ||
const size_type< matrix_t > & | n1, | ||
const size_type< matrix_t > & | n2 | ||
) |
schur_swap, swaps 2 eigenvalues of A.
[in] | want_q | bool Whether or not to apply the transformations to Q |
[in,out] | A | n-by-n matrix. Must be in Schur form |
[in,out] | Q | n-by-n matrix. Orthogonal matrix, not referenced if want_q is false |
[in] | j0 | integer Index of first eigenvalue block |
[in] | n1 | integer Size of first eigenvalue block |
[in] | n2 | integer Size of second eigenvalue block |
Implementation for complex matrices
Computes the singular value decomposition of a 2-by-2 real triangular matrix.
T = [ F G ] [ 0 H ].
On return, SSMAX is the larger singular value and SSMIN is the smaller singular value.
[in] | f | scalar, T(0,0). |
[in] | g | scalar, T(0,1). |
[in] | h | scalar, T(1,1). |
[out] | ssmin | scalar. ssmin is the smaller singular value. |
[out] | ssmax | scalar. ssmax is the larger singular value. |
void tlapack::svd22 | ( | const T & | f, |
const T & | g, | ||
const T & | h, | ||
T & | ssmin, | ||
T & | ssmax, | ||
T & | csl, | ||
T & | snl, | ||
T & | csr, | ||
T & | snr | ||
) |
Computes the singular value decomposition of a 2-by-2 real triangular matrix.
T = [ F G ] [ 0 H ].
On return, abs(SSMAX) is the larger singular value, abs(SSMIN) is the smaller singular value, and (CSL,SNL) and (CSR,SNR) are the left and right singular vectors for abs(SSMAX), giving the decomposition
[ CSL SNL ] [ F G ] [ CSR -SNR ] = [ SSMAX 0 ] [-SNL CSL ] [ 0 H ] [ SNR CSR ] [ 0 SSMIN ].
[in] | f | scalar, T(0,0). |
[in] | g | scalar, T(0,1). |
[in] | h | scalar, T(1,1). |
[out] | ssmin | scalar. abs(ssmin) is the smaller singular value. |
[out] | ssmax | scalar. abs(ssmax) is the larger singular value. |
[out] | csl | scalar. Cosine factor of the rotation from the left. |
[out] | snl | scalar. Sine factor of the rotation from the left. The vector (CSL, SNL) is a unit left singular vector for the singular value abs(SSMAX). |
[out] | csr | scalar. Cosine factor of the rotation from the right. |
[out] | snr | scalar. Sine factor of the rotation from the right. The vector (CSR, SNR) is a unit right singular vector for the singular value abs(SSMAX). |
void tlapack::transpose | ( | matrixA_t & | A, |
matrixB_t & | B, | ||
const TransposeOpts & | opts = {} |
||
) |
transpose a matrix A into a matrix B.
[in] | A | m-by-n matrix The matrix to be transposed |
[out] | B | n-by-m matrix On exit, B = A**T |
[in] | opts | Options. |
int tlapack::ul_mult | ( | matrix_t & | A | ) |
ul_mult computes the matrix product of an upper triangular matrix U and a lower triangular unital matrix L Given input matrix A, nonzero part of L is the subdiagonal of A and on the diagonal of L is assumed to be 1, and the nonzero part of U is diagonal and super-diagonal part of A
[in,out] | A | n-by-n complex matrix. On entry, subdiagonal of A contains L(lower triangular and unital) and diagonal and superdiagonal part of A contains U(upper triangular). On exit, A is overwritten by L*U |