<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
Auxiliary routines

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>
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, TYtlapack::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, TZtlapack::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
 

Detailed Description

Auxiliary routines that are used by the computational routines.


Function Documentation

◆ check_generalized_similarity_transform() [1/2]

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

Returns
frobenius norm of res
Parameters
[in]An by n matrix
[in]Qn by n unitary matrix
[in]Zn by n unitary matrix
[in]Bn by n matrix

◆ check_generalized_similarity_transform() [2/2]

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.

Returns
frobenius norm of res
Parameters
[in]An by n matrix
[in]Qn by n unitary matrix
[in]Zn by n unitary matrix
[in]Bn by n matrix
[out]resn by n matrix as defined above
[out]workn by n workspace matrix

◆ check_orthogonality() [1/2]

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.

Returns
frobenius norm of error
Parameters
[in]Qm by n (almost) orthogonal matrix

◆ check_orthogonality() [2/2]

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.

Returns
frobenius norm of res
Parameters
[in]Qm by n (almost) orthogonal matrix
[out]resn by n matrix as defined above

◆ check_similarity_transform() [1/2]

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

Returns
frobenius norm of res
Parameters
[in]An by n matrix
[in]Qn by n unitary matrix
[in]Bn by n matrix

◆ check_similarity_transform() [2/2]

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.

Returns
frobenius norm of res
Parameters
[in]An by n matrix
[in]Qn by n unitary matrix
[in]Bn by n matrix
[out]resn by n matrix as defined above
[out]workn by n workspace matrix

◆ conjtranspose()

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.

Parameters
[in]Am-by-n matrix The matrix to be transposed
[out]Bn-by-m matrix On exit, B = A**H
[in]optsOptions.

◆ conjugate()

template<TLAPACK_VECTOR vector_t>
void tlapack::conjugate ( vector_t x)

Conjugates a vector.

Parameters
[in]xVector of size n.

◆ generalized_schur_move()

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.

Returns
0 if success
1 two adjacent blocks were too close to swap (the problem is very ill-conditioned); the pencil may have been partially reordered, and ILST points to the first row of the current position of the block being moved.
Parameters
[in]want_qbool Whether or not to apply the transformations to Q
[in]want_zbool Whether or not to apply the transformations to Q
[in,out]An-by-n matrix. Must be in Schur form
[in,out]Bn-by-n matrix. Must be in Schur form
[in,out]Qn-by-n matrix. unitary matrix, not referenced if want_q is false
[in,out]Zn-by-n matrix. unitary matrix, not referenced if want_q is false
[in,out]ifstinteger Initial row index of the eigenvalue block
[in,out]ilstinteger 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.

◆ generalized_schur_swap()

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.

Returns
0 if success
1 the swap failed, this usually means the eigenvalues of the blocks are too close.
Parameters
[in]want_qbool Whether or not to apply the transformations to Q
[in]want_zbool Whether or not to apply the transformations to Z
[in,out]An-by-n matrix. Must be in Schur form
[in,out]Bn-by-n matrix. Must be in Schur form
[in,out]Qn-by-n matrix. Orthogonal matrix, not referenced if want_q is false
[in,out]Zn-by-n matrix. Orthogonal matrix, not referenced if want_q is false
[in]j0integer Index of first eigenvalue block
[in]n1integer Size of first eigenvalue block
[in]n2integer Size of second eigenvalue block

Implementation for complex matrices

◆ infnorm_colmajor()

template<TLAPACK_MATRIX matrix_t>
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.

See also
lange() for the generic implementation that does not use workspaces.
Parameters
[in]Am-by-n matrix.

◆ infnorm_colmajor_work()

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.

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.

See also
lange() for the generic implementation that does not use workspaces.
Parameters
[in]Am-by-n matrix.
workWorkspace. Use the workspace query to determine the size needed.

◆ infnorm_hermitian_colmajor()

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.

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.

See also
lanhe() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]An-by-n matrix.

◆ infnorm_hermitian_colmajor_work()

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.

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.

See also
lanhe() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]An-by-n matrix.
workWorkspace. Use the workspace query to determine the size needed.

◆ infnorm_symmetric_colmajor()

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.

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.

See also
lansy() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]An-by-n matrix.

◆ infnorm_symmetric_colmajor_work()

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.

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.

See also
lansy() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]An-by-n matrix.
workWorkspace. Use the workspace query to determine the size needed.

◆ infnorm_triangular_colmajor()

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.

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.

See also
lantr() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]diagWhether A has a unit or non-unit diagonal:
[in]Am-by-n matrix.

◆ infnorm_triangular_colmajor_work()

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.

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.

See also
lantr() for the generic implementation that does not use workspaces.
Parameters
[in]uplo
[in]diagWhether A has a unit or non-unit diagonal:
[in]Am-by-n matrix.
workWorkspace. Use the workspace query to determine the size needed.

◆ inv_house3()

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.

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

Parameters
[in]A2x3 matrix.
[out]vvector of size 3
[out]tau

◆ labrd()

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.

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

Returns
0 if success
Parameters
[in,out]Am-by-n matrix.
[out]tauqvector of length nb. The scalar factors of the elementary reflectors which represent the unitary matrix Q.
[out]taupvector of length nb. The scalar factors of the elementary reflectors which represent the unitary matrix P.
[out]Xm-by-nb matrix.
[out]Yn-by-nb matrix.

◆ lacpy()

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 Parameters
uplo_tEither Uplo or any class that implements operator Uplo().
Parameters
[in]uplo
  • Uplo::Upper: Upper triangle of A and B are referenced;
  • Uplo::Lower: Lower triangle of A and B are referenced;
  • Uplo::General: All entries of A are referenced; the first m rows of B and first n columns of B are referenced.
Am-by-n matrix.
Bmatrix with at least m rows and at least n columns.

◆ ladiv() [1/2]

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.

\[ p + iq = (a + ib) / (c + id) \]

Template Parameters
real_tFloating-point type.
Parameters
[in]aReal part of numerator.
[in]bImaginary part of numerator.
[in]cReal part of denominator.
[in]dImaginary part of denominator.
[out]pReal part of quotient.
[out]qImaginary part of quotient.

◆ ladiv() [2/2]

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.

Returns
x/y
Template Parameters
real_tFloating-point type.
Parameters
[in]xComplex numerator.
[in]yComplex denominator.

◆ lahqr()

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.

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

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

Implementation for complex matrices.

◆ lahqr_eig22()

template<TLAPACK_SCALAR T>
void tlapack::lahqr_eig22 ( a00,
a01,
a10,
a11,
complex_type< T > &  s1,
complex_type< T > &  s2 
)

Computes the eigenvalues of a 2x2 matrix A.

Parameters
[in]a00Element (0,0) of A.
[in]a01Element (0,1) of A.
[in]a10Element (1,0) of A.
[in]a11Element (1,1) of A.
[out]s1
[out]s2s1 and s2 are the eigenvalues of A

◆ lahqr_schur22()

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.

A = [a b] = [cs -sn] [aa bb] [ cs sn] [c d] [sn cs] [cc dd] [-sn cs]

Parameters
[in,out]ascalar, A(0,0).
[in,out]bscalar, A(0,1).
[in,out]cscalar, A(1,0).
[in,out]dscalar, 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]s1complex scalar. First eigenvalue of the matrix.
[out]s2complex scalar. Second eigenvalue of the matrix.
[out]csscalar. Cosine factor of the rotation
[out]snscalar. Sine factor of the rotation

◆ lahqr_shiftcolumn() [1/2]

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.

This is used to introduce shifts in the QR algorithm

Returns
0 if success
Parameters
[in]H2x2 or 3x3 matrix. The matrix H as in the formula above.
[out]vvector of size 2 or 3 On exit, a multiple of the product
[in]s1The scalar s1 as in the formula above
[in]s2The scalar s2 as in the formula above

◆ lahqr_shiftcolumn() [2/2]

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.

This is used to introduce shifts in the QR algorithm

Returns
0 if success
Parameters
[in]H2x2 or 3x3 matrix. The matrix H as in the formula above.
[out]vvector of size 2 or 3 On exit, a multiple of the product
[in]s1The scalar s1 as in the formula above
[in]s2The scalar s2 as in the formula above

◆ lahqz()

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.

Returns
0 if success
i if the QZ algorithm failed to compute all the eigenvalues in a total of 30 iterations per eigenvalue. elements i:ihi of alpha and beta contain those eigenvalues which have been successfully computed.
Parameters
[in]want_sbool. If true, the full Schur form will be computed.
[in]want_qbool. If true, the Schur vectors Q will be computed.
[in]want_zbool. If true, the Schur vectors Z will be computed.
[in]ilointeger. Either ilo=0 or A(ilo,ilo-1) = 0.
[in]ihiinteger. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n.
[in,out]An by n matrix.
[in,out]Bn by n matrix.
[out]alphasize n vector.
[out]betasize n vector.
[in,out]Qn by n matrix.
[in,out]Zn by n matrix.

◆ lahqz_shiftcolumn()

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.

This is used to introduce shifts in the QZ algorithm

Returns
0 if success
Parameters
[in]A2x2 or 3x3 matrix. The matrix A as in the formula above.
[in]B2x2 or 3x3 matrix. The matrix B as in the formula above.
[out]vvector of size 2 or 3 On exit, a multiple of the product
[in]s1The scalar s1 as in the formula above
[in]s2The scalar s2 as in the formula above
[in]beta1The scalar beta1 as in the formula above
[in]beta2The scalar beta2 as in the formula above

◆ lahr2()

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.

The matrix Q is represented as a product of elementary reflectors

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

Each H_i has the form

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

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

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

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

Returns
0 if success
Parameters
[in]kinteger
[in]nbinteger
[in,out]An-by-n matrix.
[out]tauReal vector of length n-1. The scalar factors of the elementary reflectors.
[out]Tnb-by-nb matrix.
[out]Yn-by-nb matrix.

◆ lange()

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 Parameters
norm_tEither Norm or any class that implements operator Norm().
Parameters
[in]normType
  • Norm::Max: Maximum absolute value over all elements of the matrix. Note: this is not a consistent matrix norm.
  • Norm::One: 1-norm, the maximum value of the absolute sum of each column.
  • Norm::Inf: Inf-norm, the maximum value of the absolute sum of each row.
  • Norm::Fro: Frobenius norm of the matrix. Square root of the sum of the square of each entry in the matrix.
[in]Am-by-n matrix.

◆ lanhe()

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 Parameters
norm_tEither Norm or any class that implements operator Norm().
uplo_tEither Uplo or any class that implements operator Uplo().
Parameters
[in]normType
  • Norm::Max: Maximum absolute value over all elements of the matrix. Note: this is not a consistent matrix norm.
  • Norm::One: 1-norm, the maximum value of the absolute sum of each column.
  • Norm::Inf: Inf-norm, the maximum value of the absolute sum of each row.
  • Norm::Fro: Frobenius norm of the matrix. Square root of the sum of the square of each entry in the matrix.
[in]uplo
[in]An-by-n hermitian matrix.

◆ lansy()

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 Parameters
norm_tEither Norm or any class that implements operator Norm().
uplo_tEither Uplo or any class that implements operator Uplo().
Parameters
[in]normType
  • Norm::Max: Maximum absolute value over all elements of the matrix. Note: this is not a consistent matrix norm.
  • Norm::One: 1-norm, the maximum value of the absolute sum of each column.
  • Norm::Inf: Inf-norm, the maximum value of the absolute sum of each row.
  • Norm::Fro: Frobenius norm of the matrix. Square root of the sum of the square of each entry in the matrix.
[in]uplo
[in]An-by-n symmetric matrix.

◆ lantr()

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 Parameters
norm_tEither Norm or any class that implements operator Norm().
uplo_tEither Uplo or any class that implements operator Uplo().
diag_tEither Diag or any class that implements operator Diag().
Parameters
[in]normType
  • Norm::Max: Maximum absolute value over all elements of the matrix. Note: this is not a consistent matrix norm.
  • Norm::One: 1-norm, the maximum value of the absolute sum of each column.
  • Norm::Inf: Inf-norm, the maximum value of the absolute sum of each row.
  • Norm::Fro: Frobenius norm of the matrix. Square root of the sum of the square of each entry in the matrix.
[in]uplo
[in]diagWhether A has a unit or non-unit diagonal:
[in]Am-by-n triangular matrix.

◆ lapy2()

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.

Returns
\(\sqrt{x^2+y^2}\)
Parameters
[in]xscalar value x
[in]yscalar value y

◆ lapy3()

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.

Returns
\(\sqrt{x^2+y^2+z^2}\)
Template Parameters
real_tFloating-point type.
Parameters
[in]xscalar value x
[in]yscalar value y
[in]zscalar value y

◆ larf() [1/2]

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.

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.

Template Parameters
side_tEither Side or any class that implements operator Side().
Parameters
[in]side
[in]directionv = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]vVector of size m if side = Side::Left, or n if side = Side::Right.
[in]tauValue of tau in the representation of H.
[in,out]COn 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.

◆ larf() [2/2]

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.

\[ 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} \]

Template Parameters
side_tEither Side or any class that implements operator Side().
Parameters
[in]side
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]xVector \(v\) if storeMode = StoreV::Columnwise, or \(v^H\) if storeMode = StoreV::Rowwise.
[in]tauValue of tau in the representation of H.
[in,out]C0Vector of size n if side = Side::Left, or m if side = Side::Right. On exit, C0 is overwritten by
  • \((1-\tau) C_0 - \tau v^H C_1\) if side = Side::Left, or
  • \((1-\tau) C_0 -\tau C_1 v\) if side = Side::Right.
[in,out]C1Matrix 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
  • \(-\tau v C_0 + (I-\tau vv^H) C_1\) if side = Side::Left, or
  • \(-\tau C_0 v^H + C_1 (I-\tau vv^H)\) if side = Side::Right.

◆ larf_work() [1/2]

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.

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.

Template Parameters
side_tEither Side or any class that implements operator Side().
Parameters
[in]side
[in]directionv = [ 1 x ] if direction == Direction::Forward and v = [ x 1 ] if direction == Direction::Backward.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]vVector of size m if side = Side::Left, or n if side = Side::Right.
[in]tauValue of tau in the representation of H.
[in,out]COn 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.
workWorkspace. Use the workspace query to determine the size needed.

◆ larf_work() [2/2]

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.

\[ 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} \]

Template Parameters
side_tEither Side or any class that implements operator Side().
Parameters
[in]side
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]xVector \(v\) if storeMode = StoreV::Columnwise, or \(v^H\) if storeMode = StoreV::Rowwise.
[in]tauValue of tau in the representation of H.
[in,out]C0Vector of size n if side = Side::Left, or m if side = Side::Right. On exit, C0 is overwritten by
  • \((1-\tau) C_0 - \tau v^H C_1\) if side = Side::Left, or
  • \((1-\tau) C_0 -\tau C_1 v\) if side = Side::Right.
[in,out]C1Matrix 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
  • \(-\tau v C_0 + (I-\tau vv^H) C_1\) if side = Side::Left, or
  • \(-\tau C_0 v^H + C_1 (I-\tau vv^H)\) if side = Side::Right.
workWorkspace. Use the workspace query to determine the size needed.

◆ larfb()

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.

Parameters
[in]side
[in]trans
  • Op::NoTrans: apply \(H \) (No transpose).
  • Op::Trans: apply \(H^T\) (Transpose, only allowed if the type of H is Real).
  • Op::ConjTrans: apply \(H^H\) (Conjugate transpose).
[in]directionIndicates how H is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]V
[in]TmatrixThe k-by-k matrix T. The triangular k-by-k matrix T in the representation of the block reflector.
[in,out]COn 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\).
Further Details

The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with 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 )

◆ larfb_work()

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.

Parameters
[in]side
[in]trans
  • Op::NoTrans: apply \(H \) (No transpose).
  • Op::Trans: apply \(H^T\) (Transpose, only allowed if the type of H is Real).
  • Op::ConjTrans: apply \(H^H\) (Conjugate transpose).
[in]directionIndicates how H is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]V
[in]TmatrixThe k-by-k matrix T. The triangular k-by-k matrix T in the representation of the block reflector.
[in,out]COn 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\).
Further Details

The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with 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 )
Parameters
workWorkspace. Use the workspace query to determine the size needed.

◆ larfg() [1/2]

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.

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.

Template Parameters
direction_tEither Direction or any class that implements operator Direction().
storage_tEither StoreV or any class that implements operator StoreV().
Parameters
[in]directionv = [ alpha x ] if direction == Direction::Forward and v = [ x alpha ] if direction == Direction::Backward.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in,out]vVector 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]tauThe value tau.

◆ larfg() [2/2]

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.

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.

Parameters
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in,out]alphaOn entry, the value alpha. On exit, it is overwritten with the value beta.
[in,out]xVector of length n-1. On entry, the vector x. On exit, it is overwritten with the vector y.
[out]tauThe value tau.

◆ larft()

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.

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 )
Template Parameters
direction_tEither Direction or any class that implements operator Direction().
storage_tEither StoreV or any class that implements operator StoreV().
Parameters
[in]directionIndicates how H is formed from a product of elementary reflectors.
[in]storeModeIndicates how the vectors which define the elementary reflectors are stored:
[in]V
[in]tauVector of length k containing the scalar factors of the elementary reflectors H.
[out]TMatrix of size k-by-k containing the triangular factors of the block reflector.

◆ larnv()

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.

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.

Template Parameters
idistSpecifies the distribution:
  1. real and imaginary parts each uniform (0,1).
  2. real and imaginary parts each uniform (-1,1).
  3. real and imaginary parts each normal (0,1).
  4. uniformly distributed on the disc abs(z) < 1.
  5. uniformly distributed on the circle abs(z) = 1.
Parameters
[in,out]iseedSeed for the random number generator. The seed is updated inside the routine ( seed := seed + 1 ).
[out]xVector of length n.

◆ lascl() [1/2]

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.

Specific implementation for band access types.

Parameters
[in]accessTypeDetermines the entries of A that are scaled by a/b.
[in]bThe denominator of the scalar a/b.
[in]aThe numerator of the scalar a/b.
[in,out]AMatrix to be scaled by a/b.
See also
lascl(uplo_t uplo, const b_type& b, const a_type& a, matrix_t& A)

◆ lascl() [2/2]

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.

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.

Template Parameters
uplo_tType of access inside the algorithm. Either Uplo or any type that implements operator Uplo().
matrix_tMatrix type.
a_typeType of the coefficient a. a_type cannot be a complex type.
b_typeType of the coefficient b. b_type cannot be a complex type.
Parameters
[in]uploDetermines 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]bThe denominator of the scalar a/b.
[in]aThe numerator of the scalar a/b.
[in,out]AMatrix to be scaled by a/b.
Returns
0 if success..

◆ laset()

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 Parameters
uplo_tEither Uplo or any class that implements operator Uplo().
Parameters
[in]uplo
[in]alphaValue to assign to the off-diagonal elements of A.
[in]betaValue to assign to the diagonal elements of A.
[out]Am-by-n matrix.

◆ lassq() [1/2]

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.

\[ scl smsq := \sum_{i = 0}^n |x_i|^2 + scale^2 sumsq. \]

See also
lassq(const vector_t& x, real_type<type_t<vector_t>>& scale, real_type<type_t<vector_t>>& sumsq, abs_f absF).

Specific implementation using

absF = []( const T& x ) { return abs( x ); }
typename traits::real_type_traits< Types..., int >::type real_type
The common real type of the list of types.
Definition scalar_type_traits.hpp:113

where T is the type_t< vector_t >.

◆ lassq() [2/2]

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.

\[ scl smsq := \sum_{i = 0}^n |x_i|^2 + scale^2 sumsq, \]

scale and sumsq are assumed to be non-negative.

Parameters
[in]xVector of size n.
[in,out]scaleReal 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]sumsqReal 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]absFLambda function that computes the absolute value.
absF = []( const T& x ) -> real_type<T> { ... };

◆ lasy2()

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.

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.

Note
Sets scale = 1 and xnorm = 0 if N1 = N2 = 0.
Todo:
Implement for n1=2 and n2=1 or n1=1 and n2=2.

◆ lu_mult()

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.

this is the recursive variant

Parameters
[in,out]An-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]optsOptions.

◆ move_bulge()

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.

If the bulge collapses, an attempt is made to reintroduce it using shifts s1 and s2.

Parameters
[in,out]H4x3 matrix.
[in,out]vvector 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]s1complex valued shift
[in]s2complex valued shift

◆ multishift_qz()

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.

Returns
0 if success
i if the QZ algorithm failed to compute all the eigenvalues in a total of 30 iterations per eigenvalue. elements i:ihi of alpha and beta contain those eigenvalues which have been successfully computed.
Parameters
[in]want_sbool. If true, the full Schur form will be computed.
[in]want_qbool. If true, the Schur vectors Q will be computed.
[in]want_zbool. If true, the Schur vectors Z will be computed.
[in]ilointeger. Either ilo=0 or A(ilo,ilo-1) = 0.
[in]ihiinteger. The matrix A is assumed to be already quasi-triangular in rows and columns ihi:n.
[in,out]An by n matrix.
[in,out]Bn by n matrix.
[out]alphasize n vector.
[out]betasize n vector.
[in,out]Qn by n matrix.
[in,out]Zn by n matrix.
[in,out]optsOptions.

◆ rscl()

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\).

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:

  1. If the real part of alpha is zero, then scale by the reciprocal of the imaginary part of alpha.
  2. If either real or imaginary part is greater than or equal to safeMax. If so, do proper scaling using a reliable algorithm for complex division.
  3. If 1 and 2 are false, we can compute the reciprocal of real and imaginary parts of 1/alpha without NaNs. If both components are in the safe range, then we can do the reciprocal without overflow or underflow. Otherwise, we do proper scaling.
Note
For complex alpha, it is important to always scale first by the smallest factor. This avoids an early overflow that may lead to a NaN. For example, if alpha = 5.877472e-39 + 2.802597e-45 * I and v[i] = 30 + 30 * I in single precision, we want to scale v[i] first by (safeMin / alpha) and after that scale the result by safeMax. Doing it the other way around will possibly lead to a (Inf - Inf) which should be tagged as a NaN.
Parameters
[in]alphaScalar.
[in,out]xA n-element vector.

◆ schur_move()

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.

Returns
0 if success
1 two adjacent blocks were too close to swap (the problem is very ill-conditioned); T may have been partially reordered, and ILST points to the first row of the current position of the block being moved.
Parameters
[in]want_qbool Whether or not to apply the transformations to Q
[in,out]An-by-n matrix. Must be in Schur form
[in,out]Qn-by-n matrix. Orthogonal matrix, not referenced if want_q is false
[in,out]ifstinteger Initial row index of the eigenvalue block
[in,out]ilstinteger 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.

◆ schur_swap()

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.

Returns
0 if success
1 the swap failed, this usually means the eigenvalues of the blocks are too close.
Parameters
[in]want_qbool Whether or not to apply the transformations to Q
[in,out]An-by-n matrix. Must be in Schur form
[in,out]Qn-by-n matrix. Orthogonal matrix, not referenced if want_q is false
[in]j0integer Index of first eigenvalue block
[in]n1integer Size of first eigenvalue block
[in]n2integer Size of second eigenvalue block

Implementation for complex matrices

◆ singularvalues22()

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.

T = [  F   G  ]
    [  0   H  ].

On return, SSMAX is the larger singular value and SSMIN is the smaller singular value.

Parameters
[in]fscalar, T(0,0).
[in]gscalar, T(0,1).
[in]hscalar, T(1,1).
[out]ssminscalar. ssmin is the smaller singular value.
[out]ssmaxscalar. ssmax is the larger singular value.

◆ svd22()

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.

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

Parameters
[in]fscalar, T(0,0).
[in]gscalar, T(0,1).
[in]hscalar, T(1,1).
[out]ssminscalar. abs(ssmin) is the smaller singular value.
[out]ssmaxscalar. abs(ssmax) is the larger singular value.
[out]cslscalar. Cosine factor of the rotation from the left.
[out]snlscalar. 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]csrscalar. Cosine factor of the rotation from the right.
[out]snrscalar. Sine factor of the rotation from the right. The vector (CSR, SNR) is a unit right singular vector for the singular value abs(SSMAX).

◆ transpose()

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.

Parameters
[in]Am-by-n matrix The matrix to be transposed
[out]Bn-by-m matrix On exit, B = A**T
[in]optsOptions.

◆ ul_mult()

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

Returns
0 if success.
Parameters
[in,out]An-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