fortran-lapack
|
Functions/Subroutines | |
pure subroutine, public | la_waxpy (n, za, zx, incx, zy, incy) |
ZAXPY: constant times a vector plus a vector. | |
pure subroutine, public | la_wcopy (n, zx, incx, zy, incy) |
ZCOPY: copies a vector, x, to a vector, y. | |
pure complex(qp) function, public | la_wdotc (n, zx, incx, zy, incy) |
ZDOTC: forms the dot product of two complex vectors ZDOTC = X^H * Y. | |
pure complex(qp) function, public | la_wdotu (n, zx, incx, zy, incy) |
ZDOTU: forms the dot product of two complex vectors ZDOTU = X^T * Y. | |
pure subroutine, public | la_wdrot (n, zx, incx, zy, incy, c, s) |
Applies a plane rotation, where the cos and sin (c and s) are real and the vectors cx and cy are complex. jack dongarra, linpack, 3/11/78. | |
pure subroutine, public | la_wdscal (n, da, zx, incx) |
ZDSCAL: scales a vector by a constant. | |
pure subroutine, public | la_wgbmv (trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy) |
ZGBMV: performs one of the matrix-vector operations y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y, or y := alpha*A**H*x + beta*y, where alpha and beta are scalars, x and y are vectors and A is an m by n band matrix, with kl sub-diagonals and ku super-diagonals. | |
pure subroutine, public | la_wgemm (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc) |
ZGEMM: performs one of the matrix-matrix operations C := alpha*op( A )*op( B ) + beta*C, where op( X ) is one of op( X ) = X or 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. | |
pure subroutine, public | la_wgemv (trans, m, n, alpha, a, lda, x, incx, beta, y, incy) |
ZGEMV: performs one of the matrix-vector operations y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y, or y := alpha*A**H*x + beta*y, where alpha and beta are scalars, x and y are vectors and A is an m by n matrix. | |
pure subroutine, public | la_wgerc (m, n, alpha, x, incx, y, incy, a, lda) |
ZGERC: performs the rank 1 operation A := alpha*x*y**H + A, where alpha is a scalar, x is an m element vector, y is an n element vector and A is an m by n matrix. | |
pure subroutine, public | la_wgeru (m, n, alpha, x, incx, y, incy, a, lda) |
ZGERU: performs the rank 1 operation A := alpha*x*y**T + A, where alpha is a scalar, x is an m element vector, y is an n element vector and A is an m by n matrix. | |
pure subroutine, public | la_whbmv (uplo, n, k, alpha, a, lda, x, incx, beta, y, incy) |
ZHBMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian band matrix, with k super-diagonals. | |
pure subroutine, public | la_whemm (side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc) |
ZHEMM: performs one of the matrix-matrix operations C := alpha*A*B + beta*C, or C := alpha*B*A + beta*C, where alpha and beta are scalars, A is an hermitian matrix and B and C are m by n matrices. | |
pure subroutine, public | la_whemv (uplo, n, alpha, a, lda, x, incx, beta, y, incy) |
ZHEMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian matrix. | |
pure subroutine, public | la_wher (uplo, n, alpha, x, incx, a, lda) |
ZHER: performs the hermitian rank 1 operation A := alpha*x*x**H + A, where alpha is a real scalar, x is an n element vector and A is an n by n hermitian matrix. | |
pure subroutine, public | la_wher2 (uplo, n, alpha, x, incx, y, incy, a, lda) |
ZHER2: performs the hermitian rank 2 operation A := alpha*x*y**H + conjg( alpha )*y*x**H + A, where alpha is a scalar, x and y are n element vectors and A is an n by n hermitian matrix. | |
pure subroutine, public | la_wher2k (uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc) |
ZHER2K: performs one of the hermitian rank 2k operations C := alpha*A*B**H + conjg( alpha )*B*A**H + beta*C, or C := alpha*A**H*B + conjg( alpha )*B**H*A + beta*C, where alpha and beta are scalars with beta real, C is an n by n hermitian matrix and A and B are n by k matrices in the first case and k by n matrices in the second case. | |
pure subroutine, public | la_wherk (uplo, trans, n, k, alpha, a, lda, beta, c, ldc) |
ZHERK: performs one of the hermitian rank k operations 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 matrix in the first case and a k by n matrix in the second case. | |
pure subroutine, public | la_whpmv (uplo, n, alpha, ap, x, incx, beta, y, incy) |
ZHPMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian matrix, supplied in packed form. | |
pure subroutine, public | la_whpr (uplo, n, alpha, x, incx, ap) |
ZHPR: performs the hermitian rank 1 operation A := alpha*x*x**H + A, where alpha is a real scalar, x is an n element vector and A is an n by n hermitian matrix, supplied in packed form. | |
pure subroutine, public | la_whpr2 (uplo, n, alpha, x, incx, y, incy, ap) |
ZHPR2: performs the hermitian rank 2 operation A := alpha*x*y**H + conjg( alpha )*y*x**H + A, where alpha is a scalar, x and y are n element vectors and A is an n by n hermitian matrix, supplied in packed form. | |
pure subroutine, public | la_wrotg (a, b, c, s) |
! | |
pure subroutine, public | la_wscal (n, za, zx, incx) |
ZSCAL: scales a vector by a constant. | |
pure subroutine, public | la_wswap (n, zx, incx, zy, incy) |
ZSWAP: interchanges two vectors. | |
pure subroutine, public | la_wsymm (side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc) |
ZSYMM: performs one of the matrix-matrix operations C := alpha*A*B + beta*C, or C := alpha*B*A + beta*C, where alpha and beta are scalars, A is a symmetric matrix and B and C are m by n matrices. | |
pure subroutine, public | la_wsyr2k (uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc) |
ZSYR2K: performs one of the symmetric rank 2k operations 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 matrices in the first case and k by n matrices in the second case. | |
pure subroutine, public | la_wsyrk (uplo, trans, n, k, alpha, a, lda, beta, c, ldc) |
ZSYRK: performs one of the symmetric rank k operations 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 matrix in the first case and a k by n matrix in the second case. | |
pure subroutine, public | la_wtbmv (uplo, trans, diag, n, k, a, lda, x, incx) |
ZTBMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals. | |
pure subroutine, public | la_wtbsv (uplo, trans, diag, n, k, a, lda, x, incx) |
ZTBSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine. | |
pure subroutine, public | la_wtpmv (uplo, trans, diag, n, ap, x, incx) |
ZTPMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form. | |
pure subroutine, public | la_wtpsv (uplo, trans, diag, n, ap, x, incx) |
ZTPSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine. | |
pure subroutine, public | la_wtrmm (side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb) |
ZTRMM: performs one of the matrix-matrix operations B := alpha*op( A )*B, or B := alpha*B*op( A ) where alpha is a scalar, B is an m by n matrix, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T or op( A ) = A**H. | |
pure subroutine, public | la_wtrmv (uplo, trans, diag, n, a, lda, x, incx) |
ZTRMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix. | |
pure subroutine, public | la_wtrsm (side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb) |
ZTRSM: solves one of the matrix equations op( A )*X = alpha*B, or X*op( A ) = alpha*B, where alpha is a scalar, X and B are m by n matrices, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T or op( A ) = A**H. The matrix X is overwritten on B. | |
pure subroutine, public | la_wtrsv (uplo, trans, diag, n, a, lda, x, incx) |
ZTRSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element 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. | |
pure subroutine, public la_blas_w::la_waxpy | ( | integer(ilp), intent(in) | n, |
complex(qp), intent(in) | za, | ||
complex(qp), dimension(*), intent(in) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(inout) | zy, | ||
integer(ilp), intent(in) | incy ) |
ZAXPY: constant times a vector plus a vector.
pure subroutine, public la_blas_w::la_wcopy | ( | integer(ilp), intent(in) | n, |
complex(qp), dimension(*), intent(in) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(out) | zy, | ||
integer(ilp), intent(in) | incy ) |
ZCOPY: copies a vector, x, to a vector, y.
pure complex(qp) function, public la_blas_w::la_wdotc | ( | integer(ilp), intent(in) | n, |
complex(qp), dimension(*), intent(in) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | zy, | ||
integer(ilp), intent(in) | incy ) |
ZDOTC: forms the dot product of two complex vectors ZDOTC = X^H * Y.
pure complex(qp) function, public la_blas_w::la_wdotu | ( | integer(ilp), intent(in) | n, |
complex(qp), dimension(*), intent(in) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | zy, | ||
integer(ilp), intent(in) | incy ) |
ZDOTU: forms the dot product of two complex vectors ZDOTU = X^T * Y.
pure subroutine, public la_blas_w::la_wdrot | ( | integer(ilp), intent(in) | n, |
complex(qp), dimension(*), intent(inout) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(inout) | zy, | ||
integer(ilp), intent(in) | incy, | ||
real(qp), intent(in) | c, | ||
real(qp), intent(in) | s ) |
Applies a plane rotation, where the cos and sin (c and s) are real and the vectors cx and cy are complex. jack dongarra, linpack, 3/11/78.
pure subroutine, public la_blas_w::la_wdscal | ( | integer(ilp), intent(in) | n, |
real(qp), intent(in) | da, | ||
complex(qp), dimension(*), intent(inout) | zx, | ||
integer(ilp), intent(in) | incx ) |
ZDSCAL: scales a vector by a constant.
pure subroutine, public la_blas_w::la_wgbmv | ( | character, intent(in) | trans, |
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | kl, | ||
integer(ilp), intent(in) | ku, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(*), intent(inout) | y, | ||
integer(ilp), intent(in) | incy ) |
ZGBMV: performs one of the matrix-vector operations y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y, or y := alpha*A**H*x + beta*y, where alpha and beta are scalars, x and y are vectors and A is an m by n band matrix, with kl sub-diagonals and ku super-diagonals.
pure subroutine, public la_blas_w::la_wgemm | ( | character, intent(in) | transa, |
character, intent(in) | transb, | ||
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(in) | b, | ||
integer(ilp), intent(in) | ldb, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZGEMM: performs one of the matrix-matrix operations C := alpha*op( A )*op( B ) + beta*C, where op( X ) is one of op( X ) = X or 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.
pure subroutine, public la_blas_w::la_wgemv | ( | character, intent(in) | trans, |
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(*), intent(inout) | y, | ||
integer(ilp), intent(in) | incy ) |
ZGEMV: performs one of the matrix-vector operations y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y, or y := alpha*A**H*x + beta*y, where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.
pure subroutine, public la_blas_w::la_wgerc | ( | integer(ilp), intent(in) | m, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | y, | ||
integer(ilp), intent(in) | incy, | ||
complex(qp), dimension(lda,*), intent(inout) | a, | ||
integer(ilp), intent(in) | lda ) |
ZGERC: performs the rank 1 operation A := alpha*x*y**H + A, where alpha is a scalar, x is an m element vector, y is an n element vector and A is an m by n matrix.
pure subroutine, public la_blas_w::la_wgeru | ( | integer(ilp), intent(in) | m, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | y, | ||
integer(ilp), intent(in) | incy, | ||
complex(qp), dimension(lda,*), intent(inout) | a, | ||
integer(ilp), intent(in) | lda ) |
ZGERU: performs the rank 1 operation A := alpha*x*y**T + A, where alpha is a scalar, x is an m element vector, y is an n element vector and A is an m by n matrix.
pure subroutine, public la_blas_w::la_whbmv | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(*), intent(inout) | y, | ||
integer(ilp), intent(in) | incy ) |
ZHBMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian band matrix, with k super-diagonals.
pure subroutine, public la_blas_w::la_whemm | ( | character, intent(in) | side, |
character, intent(in) | uplo, | ||
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(in) | b, | ||
integer(ilp), intent(in) | ldb, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZHEMM: performs one of the matrix-matrix operations C := alpha*A*B + beta*C, or C := alpha*B*A + beta*C, where alpha and beta are scalars, A is an hermitian matrix and B and C are m by n matrices.
pure subroutine, public la_blas_w::la_whemv | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(*), intent(inout) | y, | ||
integer(ilp), intent(in) | incy ) |
ZHEMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian matrix.
pure subroutine, public la_blas_w::la_wher | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
real(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(lda,*), intent(inout) | a, | ||
integer(ilp), intent(in) | lda ) |
ZHER: performs the hermitian rank 1 operation A := alpha*x*x**H + A, where alpha is a real scalar, x is an n element vector and A is an n by n hermitian matrix.
pure subroutine, public la_blas_w::la_wher2 | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | y, | ||
integer(ilp), intent(in) | incy, | ||
complex(qp), dimension(lda,*), intent(inout) | a, | ||
integer(ilp), intent(in) | lda ) |
ZHER2: performs the hermitian rank 2 operation A := alpha*x*y**H + conjg( alpha )*y*x**H + A, where alpha is a scalar, x and y are n element vectors and A is an n by n hermitian matrix.
pure subroutine, public la_blas_w::la_wher2k | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(in) | b, | ||
integer(ilp), intent(in) | ldb, | ||
real(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZHER2K: performs one of the hermitian rank 2k operations C := alpha*A*B**H + conjg( alpha )*B*A**H + beta*C, or C := alpha*A**H*B + conjg( alpha )*B**H*A + beta*C, where alpha and beta are scalars with beta real, C is an n by n hermitian matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.
pure subroutine, public la_blas_w::la_wherk | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
real(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
real(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZHERK: performs one of the hermitian rank k operations 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 matrix in the first case and a k by n matrix in the second case.
pure subroutine, public la_blas_w::la_whpmv | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | ap, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(*), intent(inout) | y, | ||
integer(ilp), intent(in) | incy ) |
ZHPMV: performs the matrix-vector operation y := alpha*A*x + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n hermitian matrix, supplied in packed form.
pure subroutine, public la_blas_w::la_whpr | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
real(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(inout) | ap ) |
ZHPR: performs the hermitian rank 1 operation A := alpha*x*x**H + A, where alpha is a real scalar, x is an n element vector and A is an n by n hermitian matrix, supplied in packed form.
pure subroutine, public la_blas_w::la_whpr2 | ( | character, intent(in) | uplo, |
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(*), intent(in) | x, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(in) | y, | ||
integer(ilp), intent(in) | incy, | ||
complex(qp), dimension(*), intent(inout) | ap ) |
ZHPR2: performs the hermitian rank 2 operation A := alpha*x*y**H + conjg( alpha )*y*x**H + A, where alpha is a scalar, x and y are n element vectors and A is an n by n hermitian matrix, supplied in packed form.
pure subroutine, public la_blas_w::la_wrotg | ( | complex(qp), intent(inout) | a, |
complex(qp), intent(in) | b, | ||
real(qp), intent(out) | c, | ||
complex(qp), intent(out) | s ) |
!
The computation uses the formulas |x| = sqrt( Re(x)**2 + Im(x)**2 ) sgn(x) = x / |x| if x /= 0 = 1 if x = 0 c = |a| / sqrt(|a|**2 + |b|**2) s = sgn(a) * conjg(b) / sqrt(|a|**2 + |b|**2) When a and b are real and r /= 0, the formulas simplify to r = sgn(a)*sqrt(|a|**2 + |b|**2) c = a / r s = b / r the same as in DROTG when |a| > |b|. When |b| >= |a|, the sign of c and s will be different from those computed by DROTG if the signs of a and b are not the same.
pure subroutine, public la_blas_w::la_wscal | ( | integer(ilp), intent(in) | n, |
complex(qp), intent(in) | za, | ||
complex(qp), dimension(*), intent(inout) | zx, | ||
integer(ilp), intent(in) | incx ) |
ZSCAL: scales a vector by a constant.
pure subroutine, public la_blas_w::la_wswap | ( | integer(ilp), intent(in) | n, |
complex(qp), dimension(*), intent(inout) | zx, | ||
integer(ilp), intent(in) | incx, | ||
complex(qp), dimension(*), intent(inout) | zy, | ||
integer(ilp), intent(in) | incy ) |
ZSWAP: interchanges two vectors.
pure subroutine, public la_blas_w::la_wsymm | ( | character, intent(in) | side, |
character, intent(in) | uplo, | ||
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(in) | b, | ||
integer(ilp), intent(in) | ldb, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZSYMM: performs one of the matrix-matrix operations C := alpha*A*B + beta*C, or C := alpha*B*A + beta*C, where alpha and beta are scalars, A is a symmetric matrix and B and C are m by n matrices.
pure subroutine, public la_blas_w::la_wsyr2k | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(in) | b, | ||
integer(ilp), intent(in) | ldb, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZSYR2K: performs one of the symmetric rank 2k operations 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 matrices in the first case and k by n matrices in the second case.
pure subroutine, public la_blas_w::la_wsyrk | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), intent(in) | beta, | ||
complex(qp), dimension(ldc,*), intent(inout) | c, | ||
integer(ilp), intent(in) | ldc ) |
ZSYRK: performs one of the symmetric rank k operations 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 matrix in the first case and a k by n matrix in the second case.
pure subroutine, public la_blas_w::la_wtbmv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTBMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals.
pure subroutine, public la_blas_w::la_wtbsv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
integer(ilp), intent(in) | k, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTBSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
pure subroutine, public la_blas_w::la_wtpmv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), dimension(*), intent(in) | ap, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTPMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form.
pure subroutine, public la_blas_w::la_wtpsv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), dimension(*), intent(in) | ap, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTPSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
pure subroutine, public la_blas_w::la_wtrmm | ( | character, intent(in) | side, |
character, intent(in) | uplo, | ||
character, intent(in) | transa, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(inout) | b, | ||
integer(ilp), intent(in) | ldb ) |
ZTRMM: performs one of the matrix-matrix operations B := alpha*op( A )*B, or B := alpha*B*op( A ) where alpha is a scalar, B is an m by n matrix, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T or op( A ) = A**H.
pure subroutine, public la_blas_w::la_wtrmv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTRMV: performs one of the matrix-vector operations x := A*x, or x := A**T*x, or x := A**H*x, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix.
pure subroutine, public la_blas_w::la_wtrsm | ( | character, intent(in) | side, |
character, intent(in) | uplo, | ||
character, intent(in) | transa, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | m, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), intent(in) | alpha, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(ldb,*), intent(inout) | b, | ||
integer(ilp), intent(in) | ldb ) |
ZTRSM: solves one of the matrix equations op( A )*X = alpha*B, or X*op( A ) = alpha*B, where alpha is a scalar, X and B are m by n matrices, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T or op( A ) = A**H. The matrix X is overwritten on B.
pure subroutine, public la_blas_w::la_wtrsv | ( | character, intent(in) | uplo, |
character, intent(in) | trans, | ||
character, intent(in) | diag, | ||
integer(ilp), intent(in) | n, | ||
complex(qp), dimension(lda,*), intent(in) | a, | ||
integer(ilp), intent(in) | lda, | ||
complex(qp), dimension(*), intent(inout) | x, | ||
integer(ilp), intent(in) | incx ) |
ZTRSV: solves one of the systems of equations A*x = b, or A**T*x = b, or A**H*x = b, where b and x are n element 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.