<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
Legacy BLAS API.
Collaboration diagram for Legacy BLAS API.:

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

Detailed Description

Function Documentation

◆ asum()

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

Returns
1-norm of vector, \(|| Re(x) ||_1 + || Im(x) ||_1 = \sum_{i=0}^{n-1} |Re(x_i)| + |Im(x_i)|\).
Parameters
[in]nNumber of elements in x. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*incx + 1.
[in]incxStride between elements of x. incx > 0.

◆ axpy()

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

Wrapper to axpy( const alpha_t& alpha, const vectorX_t& x, vectorY_t& y ).

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in]alphaScalar alpha. If alpha is zero, y is not updated.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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()

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

Wrapper to copy( const vectorX_t& x, vectorY_t& y ).

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ dot()

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 
)
Returns
dot product, \(x^H y\).
See also
dotu for unconjugated version, \(x^T y\).

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ dotu()

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 
)
Returns
unconjugated dot product, \(x^T y\).
See also
dot for conjugated version, \(x^H y\).

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ gemm()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]transAThe operation \(op(A)\) to be used:
[in]transBThe operation \(op(B)\) to be used:
[in]mNumber of rows of the matrix C and \(op(A)\). m >= 0.
[in]nNumber of columns of the matrix C and \(op(B)\). n >= 0.
[in]kNumber of columns of \(op(A)\) and rows of \(op(B)\). k >= 0.
[in]alphaScalar alpha. If alpha is zero, A and B are not accessed.
[in]A
  • If transA = NoTrans: the m-by-k matrix A, stored in an lda-by-k array [RowMajor: m-by-lda].
  • Otherwise: the k-by-m matrix A, stored in an lda-by-m array [RowMajor: k-by-lda].
[in]ldaLeading dimension of A.
  • If transA = NoTrans: lda >= max(1, m) [RowMajor: lda >= max(1, k)].
  • Otherwise: lda >= max(1, k) [RowMajor: lda >= max(1, m)].
[in]B
  • If transB = NoTrans: the k-by-n matrix B, stored in an ldb-by-n array [RowMajor: k-by-ldb].
  • Otherwise: the n-by-k matrix B, stored in an ldb-by-k array [RowMajor: n-by-ldb].
[in]ldbLeading dimension of B.
  • If transB = NoTrans: ldb >= max(1, k) [RowMajor: ldb >= max(1, n)].
  • Otherwise: ldb >= max(1, n) [RowMajor: ldb >= max(1, k)].
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe m-by-n matrix C, stored in an ldc-by-n array [RowMajor: m-by-ldc].
[in]ldcLeading dimension of C. ldc >= max(1, m) [RowMajor: ldc >= max(1, n)].

◆ gemv()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]transThe operation to be performed:
[in]mNumber of rows of the matrix A. m >= 0.
[in]nNumber of columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A and x are not accessed.
[in]AThe m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)].
[in]x
  • If trans = NoTrans: the n-element vector x, in an array of length (n-1)*abs(incx) + 1.
  • Otherwise: the m-element vector x, in an array of length (m-1)*abs(incx) + 1.
[in]incxStride 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]betaScalar beta. If beta is zero, y need not be set on input.
[in,out]y
  • If trans = NoTrans: the m-element vector y, in an array of length (m-1)*abs(incy) + 1.
  • Otherwise: the n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ ger()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]mNumber of rows of the matrix A. m >= 0.
[in]nNumber of columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe m-element vector x, in an array of length (m-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]AThe m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)].

◆ geru()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]mNumber of rows of the matrix A. m >= 0.
[in]nNumber of columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe m-element vector x, in an array of length (m-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]AThe m-by-n matrix A, stored in an lda-by-n array [RowMajor: m-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, m) [RowMajor: lda >= max(1, n)].

◆ hemm()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]sideThe side the matrix A appears on:
[in]uploWhat part of the matrix A is referenced:
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]mNumber of rows of the matrices B and C.
[in]nNumber of columns of the matrices B and C.
[in]alphaScalar alpha. If alpha is zero, A and B are not accessed.
[in]A
  • If side = Left: the m-by-m matrix A, stored in an lda-by-m array [RowMajor: m-by-lda].
  • If side = Right: the n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A.
  • If side = Left: lda >= max(1, m).
  • If side = Right: lda >= max(1, n).
[in]BThe m-by-n matrix B, stored in ldb-by-n array [RowMajor: m-by-ldb].
[in]ldbLeading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)].
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe m-by-n matrix C, stored in ldc-by-n array [RowMajor: m-by-ldc].
[in]ldcLeading dimension of C. ldc >= max(1, m) [RowMajor: ldc >= max(1, n)].

◆ hemv()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A and x are not accessed.
[in]AThe 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]ldaLeading dimension of A. lda >= max(1, n).
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]betaScalar beta. If beta is zero, y need not be set on input.
[in,out]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ her()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]AThe 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]ldaLeading dimension of A. lda >= max(1, n).

◆ her2()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]AThe 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]ldaLeading dimension of A. lda >= max(1, n).

◆ her2k()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
  • Uplo::Lower: only the lower triangular part of C is referenced.
  • Uplo::Upper: only the upper triangular part of C is referenced.
[in]transThe operation to be performed:
[in]nNumber of rows and columns of the matrix C. n >= 0.
[in]k
  • If trans = NoTrans: number of columns of the matrix A. k >= 0.
  • Otherwise: number of rows of the matrix A. k >= 0.
[in]alphaScalar alpha. If alpha is zero, A and B are not accessed.
[in]A
  • If trans = NoTrans: the n-by-k matrix A, stored in an lda-by-k array [RowMajor: n-by-lda].
  • Otherwise: the k-by-n matrix A, stored in an lda-by-n array [RowMajor: k-by-lda].
[in]ldaLeading dimension of A.
  • If trans = NoTrans: lda >= max(1, n) [RowMajor: lda >= max(1, k)],
  • Otherwise: lda >= max(1, k) [RowMajor: lda >= max(1, n)].
[in]B
  • If trans = NoTrans: the n-by-k matrix B, stored in an ldb-by-k array [RowMajor: n-by-ldb].
  • Otherwise: the k-by-n matrix B, stored in an ldb-by-n array [RowMajor: k-by-ldb].
[in]ldbLeading dimension of B.
  • If trans = NoTrans: ldb >= max(1, n) [RowMajor: ldb >= max(1, k)],
  • Otherwise: ldb >= max(1, k) [RowMajor: ldb >= max(1, n)].
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe n-by-n Hermitian matrix C, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldcLeading dimension of C. ldc >= max(1, n).

◆ herk()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
  • Uplo::Lower: only the lower triangular part of C is referenced.
  • Uplo::Upper: only the upper triangular part of C is referenced.
[in]transThe operation to be performed:
[in]nNumber of rows and columns of the matrix C. n >= 0.
[in]k
  • If trans = NoTrans: number of columns of the matrix A. k >= 0.
  • Otherwise: number of rows of the matrix A. k >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not accessed.
[in]A
  • If trans = NoTrans: the n-by-k matrix A, stored in an lda-by-k array [RowMajor: n-by-lda].
  • Otherwise: the k-by-n matrix A, stored in an lda-by-n array [RowMajor: k-by-lda].
[in]ldaLeading 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]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe n-by-n Hermitian matrix C, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldcLeading dimension of C. ldc >= max(1, n).

◆ iamax()

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

Parameters
[in]nNumber of elements in x.
[in]xThe n-element vector x, in an array of length (n-1)*incx + 1.
[in]incxStride between elements of x. incx > 0.
Returns
In priority order:
  1. 0 if n <= 0,
  2. the index of the first NAN in \(x\) if it exists and if checkNAN == true,
  3. the index of the first Infinity in \(x\) if it exists,
  4. the Index of the infinity-norm of \(x\), \(|| x ||_{inf}\), \(\arg\max_{i=0}^{n-1} \left(|Re(x_i)| + |Im(x_i)|\right)\).
See also
iamax( const vector_t& x, const EcOpts& opts = {} )

◆ nrm2()

template<typename T >
real_type< T > tlapack::legacy::nrm2 ( idx_t  n,
T const x,
int_t  incx 
)
Returns
2-norm of vector, \(|| x ||_2 = (\sum_{i=0}^{n-1} |x_i|^2)^{1/2}\).

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*incx + 1.
[in]incxStride between elements of x. incx > 0.

◆ rot()

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:

\[ \begin{bmatrix} x^T \\ y^T \end{bmatrix} = \begin{bmatrix} c & s \\ -s & c \end{bmatrix} \begin{bmatrix} x^T \\ y^T \end{bmatrix}. \]

See also
rotg to generate the rotation.

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in,out]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]cCosine of rotation; real.
[in]sSine of rotation; complex.

◆ rotg()

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:

\[ \begin{bmatrix} r \\ 0 \end{bmatrix} = \begin{bmatrix} c & s \\ -s & c \end{bmatrix} \begin{bmatrix} a \\ b \end{bmatrix}. \]

See also
rot to apply the rotation.

Generic implementation for arbitrary data types.

Parameters
[in,out]aOn entry, scalar a. On exit, set to r.
[in,out]bOn entry, scalar b. On exit, set to s, 1/c, or 0.
[out]cCosine of rotation; real.
[out]sSine 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

◆ rotm()

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:

\[ \begin{bmatrix} x^T \\ y^T \end{bmatrix} = H \begin{bmatrix} x^T \\ y^T \end{bmatrix}. \]

See also
rotmg to generate the rotation, and for fuller description.

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in,out]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]paramArray of length 5 giving parameters of modified plane rotation.

◆ rotmg()

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.

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

See also
rotm to apply the rotation.

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:

  • For flag = -1,

    \[ H = \begin{bmatrix} h_{11} & h_{12} \\ h_{21} & h_{22} \end{bmatrix} \]

  • For flag = 0,

    \[ H = \begin{bmatrix} 1 & h_{12} \\ h_{21} & 1 \end{bmatrix} \]

  • For flag = 1,

    \[ H = \begin{bmatrix} h_{11} & 1 \\ -1 & h_{22} \end{bmatrix} \]

  • For flag = -2,

    \[ H = \begin{bmatrix} 1 & 0 \\ 0 & 1 \end{bmatrix} \]

Generic implementation for arbitrary data types.

Parameters
[in,out]d1sqrt(d1) is scaling factor for vector x.
[in,out]d2sqrt(d2) is scaling factor for vector y.
[in,out]aOn entry, scalar a. On exit, set to z.
[in]bOn entry, scalar b.
[out]paramArray 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.)

◆ scal()

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

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x. n >= 0.
[in]alphaScalar alpha.
[in]xThe n-element vector x, in an array of length (n-1)*incx + 1.
[in]incxStride between elements of x. incx > 0.

◆ swap()

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

Generic implementation for arbitrary data types.

Parameters
[in]nNumber of elements in x and y. n >= 0.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ symm()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]sideThe side the matrix A appears on:
[in]uploWhat part of the matrix A is referenced:
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]mNumber of rows of the matrices B and C.
[in]nNumber of columns of the matrices B and C.
[in]alphaScalar alpha. If alpha is zero, A and B are not accessed.
[in]A
  • If side = Left: The m-by-m matrix A, stored in an lda-by-m array.
  • If side = Right: The n-by-n matrix A, stored in an lda-by-n array.
[in]ldaLeading dimension of A.
  • If side = Left: lda >= max(1, m).
  • If side = Right: lda >= max(1, n).
[in]BThe m-by-n matrix B, stored in an ldb-by-n array.
[in]ldbLeading dimension of B. ldb >= max(1, n).
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe m-by-n matrix C, stored in an lda-by-n array.
[in]ldcLeading dimension of C. ldc >= max(1, n).

◆ symv()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A and x are not accessed.
[in]AThe n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, n).
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]betaScalar beta. If beta is zero, y need not be set on input.
[in,out]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride between elements of y. incy must not be zero. If incy < 0, uses elements of y in reverse order: y(n-1), ..., y(0).

◆ syr()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]AThe n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, n).

◆ syr2()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed from symmetry.
  • Uplo::Lower: only the lower triangular part of A is referenced.
  • Uplo::Upper: only the upper triangular part of A is referenced.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not updated.
[in]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride 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]yThe n-element vector y, in an array of length (n-1)*abs(incy) + 1.
[in]incyStride 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]AThe n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, n).

◆ syr2k()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
  • Uplo::Lower: only the lower triangular part of C is referenced.
  • Uplo::Upper: only the upper triangular part of C is referenced.
[in]transThe operation to be performed:
[in]nNumber of rows and columns of the matrix C. n >= 0.
[in]k
  • If trans = NoTrans: number of columns of the matrix A. k >= 0.
  • Otherwise: number of rows of the matrix A. k >= 0.
[in]alphaScalar alpha. If alpha is zero, A and B are not accessed.
[in]A
  • If trans = NoTrans: the n-by-k matrix A, stored in an lda-by-k array [RowMajor: n-by-lda].
  • Otherwise: the k-by-n matrix A, stored in an lda-by-n array [RowMajor: k-by-lda].
[in]ldaLeading dimension of A.
  • If trans = NoTrans: lda >= max(1, n) [RowMajor: lda >= max(1, k)],
  • Otherwise: lda >= max(1, k) [RowMajor: lda >= max(1, n)].
[in]B
  • If trans = NoTrans: the n-by-k matrix B, stored in an ldb-by-k array [RowMajor: n-by-ldb].
  • Otherwise: the k-by-n matrix B, stored in an ldb-by-n array [RowMajor: k-by-ldb].
[in]ldbLeading dimension of B.
  • If trans = NoTrans: ldb >= max(1, n) [RowMajor: ldb >= max(1, k)],
  • Otherwise: ldb >= max(1, k) [RowMajor: ldb >= max(1, n)].
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe n-by-n symmetric matrix C, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldcLeading dimension of C. ldc >= max(1, n).

◆ syrk()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix C is referenced, the opposite triangle being assumed from symmetry:
  • Uplo::Lower: only the lower triangular part of C is referenced.
  • Uplo::Upper: only the upper triangular part of C is referenced.
[in]transThe operation to be performed:
[in]nNumber of rows and columns of the matrix C. n >= 0.
[in]k
  • If trans = NoTrans: number of columns of the matrix A. k >= 0.
  • Otherwise: number of rows of the matrix A. k >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not accessed.
[in]A
  • If trans = NoTrans: the n-by-k matrix A, stored in an lda-by-k array [RowMajor: n-by-lda].
  • Otherwise: the k-by-n matrix A, stored in an lda-by-n array [RowMajor: k-by-lda].
[in]ldaLeading dimension of A.
  • If trans = NoTrans: lda >= max(1, n) [RowMajor: lda >= max(1, k)],
  • Otherwise: lda >= max(1, k) [RowMajor: lda >= max(1, n)].
[in]betaScalar beta. If beta is zero, C need not be set on input.
[in]CThe n-by-n symmetric matrix C, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldcLeading dimension of C. ldc >= max(1, n).

◆ trmm()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]sideWhether \(op(A)\) is on the left or right of B:
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed to be zero:
[in]transThe form of \(op(A)\):
[in]diagWhether A has a unit or non-unit diagonal:
[in]mNumber of rows of matrix B. m >= 0.
[in]nNumber of columns of matrix B. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not accessed.
[in]A
  • If side = Left: the m-by-m matrix A, stored in an lda-by-m array [RowMajor: m-by-lda].
  • If side = Right: the n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A.
  • If side = left: lda >= max(1, m).
  • If side = right: lda >= max(1, n).
[in,out]BThe m-by-n matrix B, stored in an ldb-by-n array [RowMajor: m-by-ldb].
[in]ldbLeading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)].

◆ trmv()

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:

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

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed to be zero.
[in]transThe operation to be performed:
[in]diagWhether A has a unit or non-unit diagonal:
  • Diag::Unit: A is assumed to be unit triangular. The diagonal elements of A are not referenced.
  • Diag::NonUnit: A is not assumed to be unit triangular.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]AThe n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, n).
[in,out]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0).

◆ trsm()

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.

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

See also
latrs for a more numerically robust implementation.

Generic implementation for arbitrary data types.

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]sideWhether \(op(A)\) is on the left or right of X:
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed to be zero:
[in]transThe form of \(op(A)\):
[in]diagWhether A has a unit or non-unit diagonal:
[in]mNumber of rows of matrices B and X. m >= 0.
[in]nNumber of columns of matrices B and X. n >= 0.
[in]alphaScalar alpha. If alpha is zero, A is not accessed.
[in]A
  • If side = Left: the m-by-m matrix A, stored in an lda-by-m array [RowMajor: m-by-lda].
  • If side = Right: the n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A.
  • If side = left: lda >= max(1, m).
  • If side = right: lda >= max(1, n).
[in,out]BOn 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]ldbLeading dimension of B. ldb >= max(1, m) [RowMajor: ldb >= max(1, n)].

◆ trsv()

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.

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

See also
LAPACK's latrs for a more numerically robust implementation.

Generic implementation for arbitrary data types.

Parameters
[in]layoutMatrix storage, Layout::ColMajor or Layout::RowMajor.
[in]uploWhat part of the matrix A is referenced, the opposite triangle being assumed to be zero.
[in]transThe equation to be solved:
[in]diagWhether A has a unit or non-unit diagonal:
  • Diag::Unit: A is assumed to be unit triangular. The diagonal elements of A are not referenced.
  • Diag::NonUnit: A is not assumed to be unit triangular.
[in]nNumber of rows and columns of the matrix A. n >= 0.
[in]AThe n-by-n matrix A, stored in an lda-by-n array [RowMajor: n-by-lda].
[in]ldaLeading dimension of A. lda >= max(1, n).
[in,out]xThe n-element vector x, in an array of length (n-1)*abs(incx) + 1.
[in]incxStride between elements of x. incx must not be zero. If incx < 0, uses elements of x in reverse order: x(n-1), ..., x(0).