<T>LAPACK 0.1.2
C++ Template Linear Algebra PACKage
|
Functions | |
template<typename T > | |
real_type< T > | tlapack::legacy::asum (idx_t n, T const *x, int_t incx) |
Wrapper to asum( vector_t const& x ). | |
template<typename TX , typename TY > | |
void | tlapack::legacy::axpy (idx_t n, scalar_type< TX, TY > alpha, TX const *x, int_t incx, TY *y, int_t incy) |
Add scaled vector, \(y = \alpha x + y\). | |
template<typename TX , typename TY > | |
void | tlapack::legacy::copy (idx_t n, TX const *x, int_t incx, TY *y, int_t incy) |
Copy vector, \(y = x\). | |
template<typename TX , typename TY > | |
scalar_type< TX, TY > | tlapack::legacy::dot (idx_t n, TX const *x, int_t incx, TY const *y, int_t incy) |
template<typename TX , typename TY > | |
scalar_type< TX, TY > | tlapack::legacy::dotu (idx_t n, TX const *x, int_t incx, TY const *y, int_t incy) |
template<typename TA , typename TB , typename TC > | |
void | tlapack::legacy::gemm (Layout layout, Op transA, Op transB, idx_t m, idx_t n, idx_t k, scalar_type< TA, TB, TC > alpha, TA const *A, idx_t lda, TB const *B, idx_t ldb, scalar_type< TA, TB, TC > beta, TC *C, idx_t ldc) |
General matrix-matrix multiply: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::gemv (Layout layout, Op trans, idx_t m, idx_t n, scalar_type< TA, TX, TY > alpha, TA const *A, idx_t lda, TX const *x, int_t incx, scalar_type< TA, TX, TY > beta, TY *y, int_t incy) |
General matrix-vector multiply: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::ger (Layout layout, idx_t m, idx_t n, scalar_type< TA, TX, TY > alpha, TX const *x, int_t incx, TY const *y, int_t incy, TA *A, idx_t lda) |
General matrix rank-1 update: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::geru (Layout layout, idx_t m, idx_t n, scalar_type< TA, TX, TY > alpha, TX const *x, int_t incx, TY const *y, int_t incy, TA *A, idx_t lda) |
General matrix rank-1 update: | |
template<typename TA , typename TB , typename TC > | |
void | tlapack::legacy::hemm (Layout layout, Side side, Uplo uplo, idx_t m, idx_t n, scalar_type< TA, TB, TC > alpha, TA const *A, idx_t lda, TB const *B, idx_t ldb, scalar_type< TA, TB, TC > beta, TC *C, idx_t ldc) |
Hermitian matrix-matrix multiply: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::hemv (Layout layout, Uplo uplo, idx_t n, scalar_type< TA, TX, TY > alpha, TA const *A, idx_t lda, TX const *x, int_t incx, scalar_type< TA, TX, TY > beta, TY *y, int_t incy) |
Hermitian matrix-vector multiply: | |
template<typename TA , typename TX > | |
void | tlapack::legacy::her (Layout layout, Uplo uplo, idx_t n, real_type< TA, TX > alpha, TX const *x, int_t incx, TA *A, idx_t lda) |
Hermitian matrix rank-1 update: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::her2 (Layout layout, Uplo uplo, idx_t n, scalar_type< TA, TX, TY > alpha, TX const *x, int_t incx, TY const *y, int_t incy, TA *A, idx_t lda) |
Hermitian matrix rank-2 update: | |
template<typename TA , typename TB , typename TC > | |
void | tlapack::legacy::her2k (Layout layout, Uplo uplo, Op trans, idx_t n, idx_t k, scalar_type< TA, TB, TC > alpha, TA const *A, idx_t lda, TB const *B, idx_t ldb, real_type< TA, TB, TC > beta, TC *C, idx_t ldc) |
Hermitian rank-k update: | |
template<typename TA , typename TC > | |
void | tlapack::legacy::herk (Layout layout, Uplo uplo, Op trans, idx_t n, idx_t k, real_type< TA, TC > alpha, TA const *A, idx_t lda, real_type< TA, TC > beta, TC *C, idx_t ldc) |
Hermitian rank-k update: | |
template<typename T > | |
idx_t | tlapack::legacy::iamax (idx_t n, T const *x, int_t incx) |
Return \(\arg\max_{i=0}^{n-1} \left(|Re(x_i)| + |Im(x_i)|\right)\). | |
template<typename T > | |
real_type< T > | tlapack::legacy::nrm2 (idx_t n, T const *x, int_t incx) |
template<typename TX , typename TY > | |
void | tlapack::legacy::rot (idx_t n, TX *x, int_t incx, TY *y, int_t incy, const real_type< TX, TY > &c, const scalar_type< TX, TY > &s) |
Apply plane rotation: | |
template<typename T > | |
void | tlapack::legacy::rotg (T *a, T *b, real_type< T > *c, T *s) |
Construct plane rotation that eliminates b, such that: | |
template<typename TX , typename TY > | |
void | tlapack::legacy::rotm (idx_t n, TX *x, int_t incx, TY *y, int_t incy, scalar_type< TX, TY > const param[5]) |
Apply modified (fast) plane rotation, H: | |
template<typename real_t > | |
void | tlapack::legacy::rotmg (real_t *d1, real_t *d2, real_t *a, real_t b, real_t param[5]) |
Construct modified (fast) plane rotation, H, that eliminates b, such that. | |
template<typename TA , typename TX > | |
void | tlapack::legacy::scal (idx_t n, const TA &alpha, TX *x, int_t incx) |
Scale vector by constant, \(x = \alpha x\). | |
template<typename TX , typename TY > | |
void | tlapack::legacy::swap (idx_t n, TX *x, int_t incx, TY *y, int_t incy) |
Swap vectors, \(x <=> y\). | |
template<typename TA , typename TB , typename TC > | |
void | tlapack::legacy::symm (Layout layout, Side side, Uplo uplo, idx_t m, idx_t n, scalar_type< TA, TB, TC > alpha, TA const *A, idx_t lda, TB const *B, idx_t ldb, scalar_type< TA, TB, TC > beta, TC *C, idx_t ldc) |
Symmetric matrix-matrix multiply: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::symv (Layout layout, Uplo uplo, idx_t n, scalar_type< TA, TX, TY > alpha, TA const *A, idx_t lda, TX const *x, int_t incx, scalar_type< TA, TX, TY > beta, TY *y, int_t incy) |
Symmetric matrix-vector multiply: | |
template<typename TA , typename TX > | |
void | tlapack::legacy::syr (Layout layout, Uplo uplo, idx_t n, scalar_type< TA, TX > alpha, TX const *x, int_t incx, TA *A, idx_t lda) |
Symmetric matrix rank-1 update: | |
template<typename TA , typename TX , typename TY > | |
void | tlapack::legacy::syr2 (Layout layout, Uplo uplo, idx_t n, scalar_type< TA, TX, TY > alpha, TX const *x, int_t incx, TY const *y, int_t incy, TA *A, idx_t lda) |
Symmetric matrix rank-2 update: | |
template<typename TA , typename TB , typename TC > | |
void | tlapack::legacy::syr2k (Layout layout, Uplo uplo, Op trans, idx_t n, idx_t k, scalar_type< TA, TB, TC > alpha, TA const *A, idx_t lda, TB const *B, idx_t ldb, scalar_type< TA, TB, TC > beta, TC *C, idx_t ldc) |
Symmetric rank-k update: | |
template<typename TA , typename TC > | |
void | tlapack::legacy::syrk (Layout layout, Uplo uplo, Op trans, idx_t n, idx_t k, scalar_type< TA, TC > alpha, TA const *A, idx_t lda, scalar_type< TA, TC > beta, TC *C, idx_t ldc) |
Symmetric rank-k update: | |
template<typename TA , typename TB > | |
void | tlapack::legacy::trmm (Layout layout, Side side, Uplo uplo, Op trans, Diag diag, idx_t m, idx_t n, scalar_type< TA, TB > alpha, TA const *A, idx_t lda, TB *B, idx_t ldb) |
Triangular matrix-matrix multiply: | |
template<typename TA , typename TX > | |
void | tlapack::legacy::trmv (Layout layout, Uplo uplo, Op trans, Diag diag, idx_t n, TA const *A, idx_t lda, TX *x, int_t incx) |
Triangular matrix-vector multiply: | |
template<typename TA , typename TB > | |
void | tlapack::legacy::trsm (Layout layout, Side side, Uplo uplo, Op trans, Diag diag, idx_t m, idx_t n, scalar_type< TA, TB > alpha, TA const *A, idx_t lda, TB *B, idx_t ldb) |
Solve the triangular matrix-vector equation. | |
template<typename TA , typename TX > | |
void | tlapack::legacy::trsv (Layout layout, Uplo uplo, Op trans, Diag diag, idx_t n, TA const *A, idx_t lda, TX *x, int_t incx) |
Solve the triangular matrix-vector equation. | |
Wrapper to asum( vector_t const& x ).
[in] | n | Number of elements in x. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*incx + 1. |
[in] | incx | Stride between elements of x. incx > 0. |
void tlapack::legacy::axpy | ( | idx_t | n, |
scalar_type< TX, TY > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TY * | y, | ||
int_t | incy | ||
) |
Add scaled vector, \(y = \alpha x + y\).
Wrapper to axpy( const alpha_t& alpha, const vectorX_t& x, vectorY_t& y ).
[in] | n | Number of elements in x and y. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, y is not updated. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
Copy vector, \(y = x\).
Wrapper to copy( const vectorX_t& x, vectorY_t& y ).
[in] | n | Number of elements in x and y. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
scalar_type< TX, TY > tlapack::legacy::dot | ( | idx_t | n, |
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy | ||
) |
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x and y. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
scalar_type< TX, TY > tlapack::legacy::dotu | ( | idx_t | n, |
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy | ||
) |
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x and y. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
void tlapack::legacy::gemm | ( | Layout | layout, |
Op | transA, | ||
Op | transB, | ||
idx_t | m, | ||
idx_t | n, | ||
idx_t | k, | ||
scalar_type< TA, TB, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB const * | B, | ||
idx_t | ldb, | ||
scalar_type< TA, TB, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
General matrix-matrix multiply:
\[ C = \alpha op(A) \times op(B) + \beta C, \]
where \(op(X)\) is one of \(op(X) = X\), \(op(X) = X^T\), or \(op(X) = X^H\), alpha and beta are scalars, and A, B, and C are matrices, with \(op(A)\) an m-by-k matrix, \(op(B)\) a k-by-n matrix, and C an m-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | transA | The operation \(op(A)\) to be used:
|
[in] | transB | The operation \(op(B)\) to be used:
|
[in] | m | Number of rows of the matrix C and \(op(A)\). m >= 0. |
[in] | n | Number of columns of the matrix C and \(op(B)\). n >= 0. |
[in] | k | Number of columns of \(op(A)\) and rows of \(op(B)\). k >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A and B are not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | B |
|
[in] | ldb | Leading dimension of B.
|
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The m-by-n matrix C, stored in an ldc-by-n array [RowMajor: m-by-ldc]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, m) [RowMajor: ldc >= max(1, n)]. |
void tlapack::legacy::gemv | ( | Layout | layout, |
Op | trans, | ||
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TX const * | x, | ||
int_t | incx, | ||
scalar_type< TA, TX, TY > | beta, | ||
TY * | y, | ||
int_t | incy | ||
) |
General matrix-vector multiply:
\[ y = \alpha op(A) x + \beta y, \]
where \(op(A)\) is one of \(op(A) = A\), \(op(A) = A^T\), or \(op(A) = A^H\), alpha and beta are scalars, x and y are vectors, and A is an m-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | trans | The operation to be performed:
|
[in] | m | Number of rows of the matrix A. m >= 0. |
[in] | n | Number of columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A and x are not accessed. |
[in] | A | The m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)]. |
[in] | x |
|
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | beta | Scalar beta. If beta is zero, y need not be set on input. |
[in,out] | y |
|
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
void tlapack::legacy::ger | ( | Layout | layout, |
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy, | ||
TA * | A, | ||
idx_t | lda | ||
) |
General matrix rank-1 update:
\[ A = \alpha x y^H + A, \]
where alpha is a scalar, x and y are vectors, and A is an m-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | m | Number of rows of the matrix A. m >= 0. |
[in] | n | Number of columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The m-element vector x, in an array of length (m-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in,out] | A | The m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)]. |
void tlapack::legacy::geru | ( | Layout | layout, |
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy, | ||
TA * | A, | ||
idx_t | lda | ||
) |
General matrix rank-1 update:
\[ A = \alpha x y^T + A, \]
where alpha is a scalar, x and y are vectors, and A is an m-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | m | Number of rows of the matrix A. m >= 0. |
[in] | n | Number of columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The m-element vector x, in an array of length (m-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in,out] | A | The m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)]. |
void tlapack::legacy::hemm | ( | Layout | layout, |
Side | side, | ||
Uplo | uplo, | ||
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TB, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB const * | B, | ||
idx_t | ldb, | ||
scalar_type< TA, TB, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Hermitian matrix-matrix multiply:
\[ C = \alpha A B + \beta C, \]
or
\[ C = \alpha B A + \beta C, \]
where alpha and beta are scalars, A is an m-by-m or n-by-n Hermitian matrix, and B and C are m-by-n matrices.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | side | The side the matrix A appears on:
|
[in] | uplo | What part of the matrix A is referenced:
|
[in] | m | Number of rows of the matrices B and C. |
[in] | n | Number of columns of the matrices B and C. |
[in] | alpha | Scalar alpha. If alpha is zero, A and B are not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | B | The m-by-n matrix B, stored in ldb-by-n array [RowMajor: m-by-ldb]. |
[in] | ldb | Leading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)]. |
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The m-by-n matrix C, stored in ldc-by-n array [RowMajor: m-by-ldc]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, m) [RowMajor: ldc >= max(1, n)]. |
void tlapack::legacy::hemv | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TX const * | x, | ||
int_t | incx, | ||
scalar_type< TA, TX, TY > | beta, | ||
TY * | y, | ||
int_t | incy | ||
) |
Hermitian matrix-vector multiply:
\[ y = \alpha A x + \beta y, \]
where alpha and beta are scalars, x and y are vectors, and A is an n-by-n Hermitian matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A and x are not accessed. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. Imaginary parts of the diagonal elements need not be set, and are assumed to be zero. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | beta | Scalar beta. If beta is zero, y need not be set on input. |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
void tlapack::legacy::her | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
real_type< TA, TX > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TA * | A, | ||
idx_t | lda | ||
) |
Hermitian matrix rank-1 update:
\[ A = \alpha x x^H + A, \]
where alpha is a scalar, x is a vector, and A is an n-by-n Hermitian matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. Imaginary parts of the diagonal elements need not be set, are assumed to be zero on entry, and are set to zero on exit. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
void tlapack::legacy::her2 | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy, | ||
TA * | A, | ||
idx_t | lda | ||
) |
Hermitian matrix rank-2 update:
\[ A = \alpha x y^H + \text{conj}(\alpha) y x^H + A, \]
where alpha is a scalar, x and y are vectors, and A is an n-by-n Hermitian matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. Imaginary parts of the diagonal elements need not be set, are assumed to be zero on entry, and are set to zero on exit. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
void tlapack::legacy::her2k | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
idx_t | n, | ||
idx_t | k, | ||
scalar_type< TA, TB, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB const * | B, | ||
idx_t | ldb, | ||
real_type< TA, TB, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Hermitian rank-k update:
\[ C = \alpha A B^H + conj(\alpha) B A^H + \beta C, \]
or
\[ C = \alpha A^H B + conj(\alpha) B^H A + \beta C, \]
where alpha and beta are scalars, C is an n-by-n Hermitian matrix, and A and B are n-by-k or k-by-n matrices.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
|
[in] | trans | The operation to be performed:
|
[in] | n | Number of rows and columns of the matrix C. n >= 0. |
[in] | k |
|
[in] | alpha | Scalar alpha. If alpha is zero, A and B are not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | B |
|
[in] | ldb | Leading dimension of B.
|
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The n-by-n Hermitian matrix C, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, n). |
void tlapack::legacy::herk | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
idx_t | n, | ||
idx_t | k, | ||
real_type< TA, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
real_type< TA, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Hermitian rank-k update:
\[ C = \alpha A A^H + \beta C, \]
or
\[ C = \alpha A^H A + \beta C, \]
where alpha and beta are real scalars, C is an n-by-n Hermitian matrix, and A is an n-by-k or k-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
|
[in] | trans | The operation to be performed:
|
[in] | n | Number of rows and columns of the matrix C. n >= 0. |
[in] | k |
|
[in] | alpha | Scalar alpha. If alpha is zero, A is not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A. If trans = NoTrans: lda >= max(1, n) [RowMajor: lda >= max(1, k)], If Otherwise: lda >= max(1, k) [RowMajor: lda >= max(1, n)]. |
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The n-by-n Hermitian matrix C, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, n). |
Return \(\arg\max_{i=0}^{n-1} \left(|Re(x_i)| + |Im(x_i)|\right)\).
[in] | n | Number of elements in x. |
[in] | x | The n-element vector x, in an array of length (n-1)*incx + 1. |
[in] | incx | Stride between elements of x. incx > 0. |
NAN
in \(x\) if it exists and if checkNAN == true
,Infinity
in \(x\) if it exists,Generic implementation for arbitrary data types.
[in] | n | Number of elements in x. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*incx + 1. |
[in] | incx | Stride between elements of x. incx > 0. |
void tlapack::legacy::rot | ( | idx_t | n, |
TX * | x, | ||
int_t | incx, | ||
TY * | y, | ||
int_t | incy, | ||
const real_type< TX, TY > & | c, | ||
const scalar_type< TX, TY > & | s | ||
) |
Apply plane rotation:
\[ \begin{bmatrix} x^T \\ y^T \end{bmatrix} = \begin{bmatrix} c & s \\ -s & c \end{bmatrix} \begin{bmatrix} x^T \\ y^T \end{bmatrix}. \]
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x and y. n >= 0. |
[in,out] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in] | c | Cosine of rotation; real. |
[in] | s | Sine of rotation; complex. |
Construct plane rotation that eliminates b, such that:
\[ \begin{bmatrix} r \\ 0 \end{bmatrix} = \begin{bmatrix} c & s \\ -s & c \end{bmatrix} \begin{bmatrix} a \\ b \end{bmatrix}. \]
Generic implementation for arbitrary data types.
[in,out] | a | On entry, scalar a. On exit, set to r. |
[in,out] | b | On entry, scalar b. On exit, set to s, 1/c, or 0. |
[out] | c | Cosine of rotation; real. |
[out] | s | Sine of rotation; scalar. |
Anderson E (2017) Algorithm 978: Safe scaling in the level 1 BLAS. ACM Trans Math Softw 44:. https://doi.org/10.1145/3061665
void tlapack::legacy::rotm | ( | idx_t | n, |
TX * | x, | ||
int_t | incx, | ||
TY * | y, | ||
int_t | incy, | ||
scalar_type< TX, TY > const | param[5] | ||
) |
Apply modified (fast) plane rotation, H:
\[ \begin{bmatrix} x^T \\ y^T \end{bmatrix} = H \begin{bmatrix} x^T \\ y^T \end{bmatrix}. \]
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x and y. n >= 0. |
[in,out] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in] | param | Array of length 5 giving parameters of modified plane rotation. |
Construct modified (fast) plane rotation, H, that eliminates b, such that.
\[ \begin{bmatrix} z \\ 0 \end{bmatrix} = H \begin{bmatrix} \sqrt{d_1} & 0 \\ 0 & \sqrt{d_2} \end{bmatrix} \begin{bmatrix} a \\ b \end{bmatrix}. \]
With modified plane rotations, vectors u and v are held in factored form as
\[ \begin{bmatrix} u^T \\ v^T \end{bmatrix} = \begin{bmatrix} \sqrt{d_1} & 0 \\ 0 & \sqrt{d_2} \end{bmatrix} \begin{bmatrix} x^T \\ y^T \end{bmatrix}. \]
Application of H to vectors x and y requires 4n flops (2n mul, 2n add) instead of 6n flops (4n mul, 2n add) as in standard plane rotations.
Let param = [ flag, \(h_{11}, h_{21}, h_{12}, h_{22}\) ]. Then H has one of the following forms:
\[ H = \begin{bmatrix} h_{11} & h_{12} \\ h_{21} & h_{22} \end{bmatrix} \]
\[ H = \begin{bmatrix} 1 & h_{12} \\ h_{21} & 1 \end{bmatrix} \]
\[ H = \begin{bmatrix} h_{11} & 1 \\ -1 & h_{22} \end{bmatrix} \]
\[ H = \begin{bmatrix} 1 & 0 \\ 0 & 1 \end{bmatrix} \]
Generic implementation for arbitrary data types.
[in,out] | d1 | sqrt(d1) is scaling factor for vector x. |
[in,out] | d2 | sqrt(d2) is scaling factor for vector y. |
[in,out] | a | On entry, scalar a. On exit, set to z. |
[in] | b | On entry, scalar b. |
[out] | param | Array of length 5 giving parameters of modified plane rotation, as described above. |
Hammarling, Sven. A note on modifications to the Givens plane rotation. IMA Journal of Applied Mathematics, 13:215-218, 1974. http://dx.doi.org/10.1093/imamat/13.2.215 (Note the notation swaps u <=> x, v <=> y, d_i -> l_i.)
Scale vector by constant, \(x = \alpha x\).
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x. n >= 0. |
[in] | alpha | Scalar alpha. |
[in] | x | The n-element vector x, in an array of length (n-1)*incx + 1. |
[in] | incx | Stride between elements of x. incx > 0. |
Swap vectors, \(x <=> y\).
Generic implementation for arbitrary data types.
[in] | n | Number of elements in x and y. n >= 0. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
void tlapack::legacy::symm | ( | Layout | layout, |
Side | side, | ||
Uplo | uplo, | ||
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TB, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB const * | B, | ||
idx_t | ldb, | ||
scalar_type< TA, TB, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Symmetric matrix-matrix multiply:
\[ C = \alpha A B + \beta C, \]
or
\[ C = \alpha B A + \beta C, \]
where alpha and beta are scalars, A is an m-by-m or n-by-n symmetric matrix, and B and C are m-by-n matrices.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | side | The side the matrix A appears on:
|
[in] | uplo | What part of the matrix A is referenced:
|
[in] | m | Number of rows of the matrices B and C. |
[in] | n | Number of columns of the matrices B and C. |
[in] | alpha | Scalar alpha. If alpha is zero, A and B are not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | B | The m-by-n matrix B, stored in an ldb-by-n array. |
[in] | ldb | Leading dimension of B. ldb >= max(1, n). |
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The m-by-n matrix C, stored in an lda-by-n array. |
[in] | ldc | Leading dimension of C. ldc >= max(1, n). |
void tlapack::legacy::symv | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TX const * | x, | ||
int_t | incx, | ||
scalar_type< TA, TX, TY > | beta, | ||
TY * | y, | ||
int_t | incy | ||
) |
Symmetric matrix-vector multiply:
\[ y = \alpha A x + \beta y, \]
where alpha and beta are scalars, x and y are vectors, and A is an n-by-n symmetric matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A and x are not accessed. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | beta | Scalar beta. If beta is zero, y need not be set on input. |
[in,out] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
void tlapack::legacy::syr | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
scalar_type< TA, TX > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TA * | A, | ||
idx_t | lda | ||
) |
Symmetric matrix rank-1 update:
\[ A = \alpha x x^T + A, \]
where alpha is a scalar, x is a vector, and A is an n-by-n symmetric matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
void tlapack::legacy::syr2 | ( | Layout | layout, |
Uplo | uplo, | ||
idx_t | n, | ||
scalar_type< TA, TX, TY > | alpha, | ||
TX const * | x, | ||
int_t | incx, | ||
TY const * | y, | ||
int_t | incy, | ||
TA * | A, | ||
idx_t | lda | ||
) |
Symmetric matrix rank-2 update:
\[ A = \alpha x y^T + \alpha y x^T + A, \]
where alpha is a scalar, x and y are vectors, and A is an n-by-n symmetric matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not updated. |
[in] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
[in] | y | The n-element vector y, in an array of length (n-1)*abs(incy) + 1. |
[in] | incy | Stride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0). |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
void tlapack::legacy::syr2k | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
idx_t | n, | ||
idx_t | k, | ||
scalar_type< TA, TB, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB const * | B, | ||
idx_t | ldb, | ||
scalar_type< TA, TB, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Symmetric rank-k update:
\[ C = \alpha A B^T + \alpha B A^T + \beta C, \]
or
\[ C = \alpha A^T B + \alpha B^T A + \beta C, \]
where alpha and beta are scalars, C is an n-by-n symmetric matrix, and A and B are n-by-k or k-by-n matrices.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
|
[in] | trans | The operation to be performed:
|
[in] | n | Number of rows and columns of the matrix C. n >= 0. |
[in] | k |
|
[in] | alpha | Scalar alpha. If alpha is zero, A and B are not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | B |
|
[in] | ldb | Leading dimension of B.
|
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The n-by-n symmetric matrix C, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, n). |
void tlapack::legacy::syrk | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
idx_t | n, | ||
idx_t | k, | ||
scalar_type< TA, TC > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
scalar_type< TA, TC > | beta, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Symmetric rank-k update:
\[ C = \alpha A A^T + \beta C, \]
or
\[ C = \alpha A^T A + \beta C, \]
where alpha and beta are scalars, C is an n-by-n symmetric matrix, and A is an n-by-k or k-by-n matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
|
[in] | trans | The operation to be performed:
|
[in] | n | Number of rows and columns of the matrix C. n >= 0. |
[in] | k |
|
[in] | alpha | Scalar alpha. If alpha is zero, A is not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in] | beta | Scalar beta. If beta is zero, C need not be set on input. |
[in] | C | The n-by-n symmetric matrix C, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | ldc | Leading dimension of C. ldc >= max(1, n). |
void tlapack::legacy::trmm | ( | Layout | layout, |
Side | side, | ||
Uplo | uplo, | ||
Op | trans, | ||
Diag | diag, | ||
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TB > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB * | B, | ||
idx_t | ldb | ||
) |
Triangular matrix-matrix multiply:
\[ B = \alpha op(A) B, \]
or
\[ B = \alpha B op(A), \]
where \(op(A)\) is one of \(op(A) = A\), \(op(A) = A^T\), or \(op(A) = A^H\), B is an m-by-n matrix, and A is an m-by-m or n-by-n, unit or non-unit, upper or lower triangular matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | side | Whether \(op(A)\) is on the left or right of B:
|
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed to be zero:
|
[in] | trans | The form of \(op(A)\):
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | m | Number of rows of matrix B. m >= 0. |
[in] | n | Number of columns of matrix B. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in,out] | B | The m-by-n matrix B, stored in an ldb-by-n array [RowMajor: m-by-ldb]. |
[in] | ldb | Leading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)]. |
void tlapack::legacy::trmv | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
Diag | diag, | ||
idx_t | n, | ||
TA const * | A, | ||
idx_t | lda, | ||
TX * | x, | ||
int_t | incx | ||
) |
Triangular matrix-vector multiply:
\[ x = op(A) x, \]
where \(op(A)\) is one of \(op(A) = A\), \(op(A) = A^T\), or \(op(A) = A^H\), x is a vector, and A is an n-by-n, unit or non-unit, upper or lower triangular matrix.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed to be zero.
|
[in] | trans | The operation to be performed:
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
[in,out] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |
void tlapack::legacy::trsm | ( | Layout | layout, |
Side | side, | ||
Uplo | uplo, | ||
Op | trans, | ||
Diag | diag, | ||
idx_t | m, | ||
idx_t | n, | ||
scalar_type< TA, TB > | alpha, | ||
TA const * | A, | ||
idx_t | lda, | ||
TB * | B, | ||
idx_t | ldb | ||
) |
Solve the triangular matrix-vector equation.
\[ op(A) X = \alpha B, \]
or
\[ X op(A) = \alpha B, \]
where \(op(A)\) is one of \(op(A) = A\), \(op(A) = A^T\), or \(op(A) = A^H\), X and B are m-by-n matrices, and A is an m-by-m or n-by-n, unit or non-unit, upper or lower triangular matrix.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | side | Whether \(op(A)\) is on the left or right of X:
|
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed to be zero:
|
[in] | trans | The form of \(op(A)\):
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | m | Number of rows of matrices B and X. m >= 0. |
[in] | n | Number of columns of matrices B and X. n >= 0. |
[in] | alpha | Scalar alpha. If alpha is zero, A is not accessed. |
[in] | A |
|
[in] | lda | Leading dimension of A.
|
[in,out] | B | On entry, the m-by-n matrix B, stored in an ldb-by-n array [RowMajor: m-by-ldb]. On exit, overwritten by the solution matrix X. |
[in] | ldb | Leading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)]. |
void tlapack::legacy::trsv | ( | Layout | layout, |
Uplo | uplo, | ||
Op | trans, | ||
Diag | diag, | ||
idx_t | n, | ||
TA const * | A, | ||
idx_t | lda, | ||
TX * | x, | ||
int_t | incx | ||
) |
Solve the triangular matrix-vector equation.
\[ op(A) x = b, \]
where \(op(A)\) is one of \(op(A) = A\), \(op(A) = A^T\), or \(op(A) = A^H\), x and b are vectors, and A is an n-by-n, unit or non-unit, upper or lower triangular matrix.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Generic implementation for arbitrary data types.
[in] | layout | Matrix storage, Layout::ColMajor or Layout::RowMajor. |
[in] | uplo | What part of the matrix A is referenced, the opposite triangle being assumed to be zero.
|
[in] | trans | The equation to be solved:
|
[in] | diag | Whether A has a unit or non-unit diagonal:
|
[in] | n | Number of rows and columns of the matrix A. n >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda]. |
[in] | lda | Leading dimension of A. lda >= max(1, n). |
[in,out] | x | The n-element vector x, in an array of length (n-1)*abs(incx) + 1. |
[in] | incx | Stride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0). |