Documentation ¶
Overview ¶
Package gonum is a Go implementation of the BLAS API. This implementation panics when the input arguments are invalid as per the standard, for example if a vector increment is zero. Note that the treatment of NaN values is not specified, and differs among the BLAS implementations. gonum.org/v1/gonum/blas/blas64 provides helpful wrapper functions to the BLAS interface. The rest of this text describes the layout of the data for the input types.
Note that in the function documentation, x[i] refers to the i^th element of the vector, which will be different from the i^th element of the slice if incX != 1.
See http://www.netlib.org/lapack/explore-html/d4/de1/_l_i_c_e_n_s_e_source.html for more license information.
Vector arguments are effectively strided slices. They have two input arguments, a number of elements, n, and an increment, incX. The increment specifies the distance between elements of the vector. The actual Go slice may be longer than necessary. The increment may be positive or negative, except in functions with only a single vector argument where the increment may only be positive. If the increment is negative, s[0] is the last element in the slice. Note that this is not the same as counting backward from the end of the slice, as len(s) may be longer than necessary. So, for example, if n = 5 and incX = 3, the elements of s are
[0 * * 1 * * 2 * * 3 * * 4 * * * ...]
where ∗ elements are never accessed. If incX = -3, the same elements are accessed, just in reverse order (4, 3, 2, 1, 0).
Dense matrices are specified by a number of rows, a number of columns, and a stride. The stride specifies the number of entries in the slice between the first element of successive rows. The stride must be at least as large as the number of columns but may be longer.
[a00 ... a0n a0* ... a1stride-1 a21 ... amn am* ... amstride-1]
Thus, dense[i*ld + j] refers to the {i, j}th element of the matrix.
Symmetric and triangular matrices (non-packed) are stored identically to Dense, except that only elements in one triangle of the matrix are accessed.
Packed symmetric and packed triangular matrices are laid out with the entries condensed such that all of the unreferenced elements are removed. So, the upper triangular matrix
[ 1 2 3 0 4 5 0 0 6 ]
and the lower-triangular matrix
[ 1 0 0 2 3 0 4 5 6 ]
will both be compacted as [1 2 3 4 5 6]. The (i, j) element of the original dense matrix can be found at element i*n - (i-1)*i/2 + j for upper triangular, and at element i * (i+1) /2 + j for lower triangular.
Banded matrices are laid out in a compact format, constructed by removing the zeros in the rows and aligning the diagonals. For example, the matrix
[ 1 2 3 0 0 0 4 5 6 7 0 0 0 8 9 10 11 0 0 0 12 13 14 15 0 0 0 16 17 18 0 0 0 0 19 20 ]
implicitly becomes (∗ entries are never accessed)
[ * 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 * 19 20 * * ]
which is given to the BLAS routine as [∗ 1 2 3 4 ...].
See http://www.crest.iu.edu/research/mtl/reference/html/banded.html for more information
Index ¶
- type Implementation
- func (Implementation) Caxpy(n int, alpha complex64, x []complex64, incX int, y []complex64, incY int)
- func (Implementation) Ccopy(n int, x []complex64, incX int, y []complex64, incY int)
- func (Implementation) Cdotc(n int, x []complex64, incX int, y []complex64, incY int) complex64
- func (Implementation) Cdotu(n int, x []complex64, incX int, y []complex64, incY int) complex64
- func (Implementation) Cgbmv(trans blas.Transpose, m, n, kL, kU int, alpha complex64, a []complex64, ...)
- func (Implementation) Cgemm(tA, tB blas.Transpose, m, n, k int, alpha complex64, a []complex64, lda int, ...)
- func (Implementation) Cgemv(trans blas.Transpose, m, n int, alpha complex64, a []complex64, lda int, ...)
- func (Implementation) Cgerc(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, ...)
- func (Implementation) Cgeru(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, ...)
- func (Implementation) Chbmv(uplo blas.Uplo, n, k int, alpha complex64, a []complex64, lda int, ...)
- func (Implementation) Chemm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex64, a []complex64, ...)
- func (Implementation) Chemv(uplo blas.Uplo, n int, alpha complex64, a []complex64, lda int, x []complex64, ...)
- func (Implementation) Cher(uplo blas.Uplo, n int, alpha float32, x []complex64, incX int, a []complex64, ...)
- func (Implementation) Cher2(uplo blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, ...)
- func (Implementation) Cher2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, ...)
- func (Implementation) Cherk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha float32, a []complex64, ...)
- func (Implementation) Chpmv(uplo blas.Uplo, n int, alpha complex64, ap []complex64, x []complex64, ...)
- func (Implementation) Chpr(uplo blas.Uplo, n int, alpha float32, x []complex64, incX int, ap []complex64)
- func (Implementation) Chpr2(uplo blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, ...)
- func (Implementation) Cscal(n int, alpha complex64, x []complex64, incX int)
- func (Implementation) Csscal(n int, alpha float32, x []complex64, incX int)
- func (Implementation) Cswap(n int, x []complex64, incX int, y []complex64, incY int)
- func (Implementation) Csymm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex64, a []complex64, ...)
- func (Implementation) Csyr2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, ...)
- func (Implementation) Csyrk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, ...)
- func (Implementation) Ctbmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex64, ...)
- func (Implementation) Ctbsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex64, ...)
- func (Implementation) Ctpmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex64, ...)
- func (Implementation) Ctpsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex64, ...)
- func (Implementation) Ctrmm(side blas.Side, uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, m, n int, ...)
- func (Implementation) Ctrmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex64, ...)
- func (Implementation) Ctrsm(side blas.Side, uplo blas.Uplo, transA blas.Transpose, diag blas.Diag, ...)
- func (Implementation) Ctrsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex64, ...)
- func (Implementation) Dasum(n int, x []float64, incX int) float64
- func (Implementation) Daxpy(n int, alpha float64, x []float64, incX int, y []float64, incY int)
- func (Implementation) Dcopy(n int, x []float64, incX int, y []float64, incY int)
- func (Implementation) Ddot(n int, x []float64, incX int, y []float64, incY int) float64
- func (Implementation) Dgbmv(tA blas.Transpose, m, n, kL, kU int, alpha float64, a []float64, lda int, ...)
- func (Implementation) Dgemm(tA, tB blas.Transpose, m, n, k int, alpha float64, a []float64, lda int, ...)
- func (Implementation) Dgemv(tA blas.Transpose, m, n int, alpha float64, a []float64, lda int, x []float64, ...)
- func (Implementation) Dger(m, n int, alpha float64, x []float64, incX int, y []float64, incY int, ...)
- func (Implementation) Dnrm2(n int, x []float64, incX int) float64
- func (Implementation) Drot(n int, x []float64, incX int, y []float64, incY int, c float64, s float64)
- func (Implementation) Drotg(a, b float64) (c, s, r, z float64)
- func (Implementation) Drotm(n int, x []float64, incX int, y []float64, incY int, p blas.DrotmParams)
- func (Implementation) Drotmg(d1, d2, x1, y1 float64) (p blas.DrotmParams, rd1, rd2, rx1 float64)
- func (Implementation) Dsbmv(ul blas.Uplo, n, k int, alpha float64, a []float64, lda int, x []float64, ...)
- func (Implementation) Dscal(n int, alpha float64, x []float64, incX int)
- func (Implementation) Dsdot(n int, x []float32, incX int, y []float32, incY int) float64
- func (Implementation) Dspmv(ul blas.Uplo, n int, alpha float64, ap []float64, x []float64, incX int, ...)
- func (Implementation) Dspr(ul blas.Uplo, n int, alpha float64, x []float64, incX int, ap []float64)
- func (Implementation) Dspr2(ul blas.Uplo, n int, alpha float64, x []float64, incX int, y []float64, ...)
- func (Implementation) Dswap(n int, x []float64, incX int, y []float64, incY int)
- func (Implementation) Dsymm(s blas.Side, ul blas.Uplo, m, n int, alpha float64, a []float64, lda int, ...)
- func (Implementation) Dsymv(ul blas.Uplo, n int, alpha float64, a []float64, lda int, x []float64, ...)
- func (Implementation) Dsyr(ul blas.Uplo, n int, alpha float64, x []float64, incX int, a []float64, ...)
- func (Implementation) Dsyr2(ul blas.Uplo, n int, alpha float64, x []float64, incX int, y []float64, ...)
- func (Implementation) Dsyr2k(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float64, a []float64, lda int, ...)
- func (Implementation) Dsyrk(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float64, a []float64, lda int, ...)
- func (Implementation) Dtbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float64, lda int, ...)
- func (Implementation) Dtbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float64, lda int, ...)
- func (Implementation) Dtpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float64, x []float64, ...)
- func (Implementation) Dtpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float64, x []float64, ...)
- func (Implementation) Dtrmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, ...)
- func (Implementation) Dtrmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float64, lda int, ...)
- func (Implementation) Dtrsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, ...)
- func (Implementation) Dtrsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float64, lda int, ...)
- func (Implementation) Dzasum(n int, x []complex128, incX int) float64
- func (Implementation) Dznrm2(n int, x []complex128, incX int) float64
- func (Implementation) Icamax(n int, x []complex64, incX int) int
- func (Implementation) Idamax(n int, x []float64, incX int) int
- func (Implementation) Isamax(n int, x []float32, incX int) int
- func (Implementation) Izamax(n int, x []complex128, incX int) int
- func (Implementation) Sasum(n int, x []float32, incX int) float32
- func (Implementation) Saxpy(n int, alpha float32, x []float32, incX int, y []float32, incY int)
- func (Implementation) Scasum(n int, x []complex64, incX int) float32
- func (Implementation) Scnrm2(n int, x []complex64, incX int) float32
- func (Implementation) Scopy(n int, x []float32, incX int, y []float32, incY int)
- func (Implementation) Sdot(n int, x []float32, incX int, y []float32, incY int) float32
- func (Implementation) Sdsdot(n int, alpha float32, x []float32, incX int, y []float32, incY int) float32
- func (Implementation) Sgbmv(tA blas.Transpose, m, n, kL, kU int, alpha float32, a []float32, lda int, ...)
- func (Implementation) Sgemm(tA, tB blas.Transpose, m, n, k int, alpha float32, a []float32, lda int, ...)
- func (Implementation) Sgemv(tA blas.Transpose, m, n int, alpha float32, a []float32, lda int, x []float32, ...)
- func (Implementation) Sger(m, n int, alpha float32, x []float32, incX int, y []float32, incY int, ...)
- func (Implementation) Snrm2(n int, x []float32, incX int) float32
- func (Implementation) Srot(n int, x []float32, incX int, y []float32, incY int, c float32, s float32)
- func (Implementation) Srotg(a, b float32) (c, s, r, z float32)
- func (Implementation) Srotm(n int, x []float32, incX int, y []float32, incY int, p blas.SrotmParams)
- func (Implementation) Srotmg(d1, d2, x1, y1 float32) (p blas.SrotmParams, rd1, rd2, rx1 float32)
- func (Implementation) Ssbmv(ul blas.Uplo, n, k int, alpha float32, a []float32, lda int, x []float32, ...)
- func (Implementation) Sscal(n int, alpha float32, x []float32, incX int)
- func (Implementation) Sspmv(ul blas.Uplo, n int, alpha float32, ap []float32, x []float32, incX int, ...)
- func (Implementation) Sspr(ul blas.Uplo, n int, alpha float32, x []float32, incX int, ap []float32)
- func (Implementation) Sspr2(ul blas.Uplo, n int, alpha float32, x []float32, incX int, y []float32, ...)
- func (Implementation) Sswap(n int, x []float32, incX int, y []float32, incY int)
- func (Implementation) Ssymm(s blas.Side, ul blas.Uplo, m, n int, alpha float32, a []float32, lda int, ...)
- func (Implementation) Ssymv(ul blas.Uplo, n int, alpha float32, a []float32, lda int, x []float32, ...)
- func (Implementation) Ssyr(ul blas.Uplo, n int, alpha float32, x []float32, incX int, a []float32, ...)
- func (Implementation) Ssyr2(ul blas.Uplo, n int, alpha float32, x []float32, incX int, y []float32, ...)
- func (Implementation) Ssyr2k(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float32, a []float32, lda int, ...)
- func (Implementation) Ssyrk(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float32, a []float32, lda int, ...)
- func (Implementation) Stbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float32, lda int, ...)
- func (Implementation) Stbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float32, lda int, ...)
- func (Implementation) Stpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float32, x []float32, ...)
- func (Implementation) Stpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float32, x []float32, ...)
- func (Implementation) Strmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, ...)
- func (Implementation) Strmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float32, lda int, ...)
- func (Implementation) Strsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, ...)
- func (Implementation) Strsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float32, lda int, ...)
- func (Implementation) Zaxpy(n int, alpha complex128, x []complex128, incX int, y []complex128, incY int)
- func (Implementation) Zcopy(n int, x []complex128, incX int, y []complex128, incY int)
- func (Implementation) Zdotc(n int, x []complex128, incX int, y []complex128, incY int) complex128
- func (Implementation) Zdotu(n int, x []complex128, incX int, y []complex128, incY int) complex128
- func (Implementation) Zdscal(n int, alpha float64, x []complex128, incX int)
- func (Implementation) Zgbmv(trans blas.Transpose, m, n, kL, kU int, alpha complex128, a []complex128, ...)
- func (Implementation) Zgemm(tA, tB blas.Transpose, m, n, k int, alpha complex128, a []complex128, lda int, ...)
- func (Implementation) Zgemv(trans blas.Transpose, m, n int, alpha complex128, a []complex128, lda int, ...)
- func (Implementation) Zgerc(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, ...)
- func (Implementation) Zgeru(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, ...)
- func (Implementation) Zhbmv(uplo blas.Uplo, n, k int, alpha complex128, a []complex128, lda int, ...)
- func (Implementation) Zhemm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex128, a []complex128, ...)
- func (Implementation) Zhemv(uplo blas.Uplo, n int, alpha complex128, a []complex128, lda int, ...)
- func (Implementation) Zher(uplo blas.Uplo, n int, alpha float64, x []complex128, incX int, a []complex128, ...)
- func (Implementation) Zher2(uplo blas.Uplo, n int, alpha complex128, x []complex128, incX int, ...)
- func (Implementation) Zher2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, ...)
- func (Implementation) Zherk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha float64, a []complex128, ...)
- func (Implementation) Zhpmv(uplo blas.Uplo, n int, alpha complex128, ap []complex128, x []complex128, ...)
- func (Implementation) Zhpr(uplo blas.Uplo, n int, alpha float64, x []complex128, incX int, ...)
- func (Implementation) Zhpr2(uplo blas.Uplo, n int, alpha complex128, x []complex128, incX int, ...)
- func (Implementation) Zscal(n int, alpha complex128, x []complex128, incX int)
- func (Implementation) Zswap(n int, x []complex128, incX int, y []complex128, incY int)
- func (Implementation) Zsymm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex128, a []complex128, ...)
- func (Implementation) Zsyr2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, ...)
- func (Implementation) Zsyrk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, ...)
- func (Implementation) Ztbmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex128, ...)
- func (Implementation) Ztbsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex128, ...)
- func (Implementation) Ztpmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex128, ...)
- func (Implementation) Ztpsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex128, ...)
- func (Implementation) Ztrmm(side blas.Side, uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, m, n int, ...)
- func (Implementation) Ztrmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex128, ...)
- func (Implementation) Ztrsm(side blas.Side, uplo blas.Uplo, transA blas.Transpose, diag blas.Diag, ...)
- func (Implementation) Ztrsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex128, ...)
Constants ¶
This section is empty.
Variables ¶
This section is empty.
Functions ¶
This section is empty.
Types ¶
type Implementation ¶
type Implementation struct{}
func (Implementation) Caxpy ¶
func (Implementation) Caxpy(n int, alpha complex64, x []complex64, incX int, y []complex64, incY int)
Caxpy adds alpha times x to y:
y[i] += alpha * x[i] for all i
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ccopy ¶
Ccopy copies the vector x to vector y.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cdotc ¶
Cdotc computes the dot product
xᴴ · y
of two complex vectors x and y.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cdotu ¶
Cdotu computes the dot product
xᵀ · y
of two complex vectors x and y.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cgbmv ¶
func (Implementation) Cgbmv(trans blas.Transpose, m, n, kL, kU int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)
Cgbmv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if trans = blas.NoTrans y = alpha * Aᵀ * x + beta * y if trans = blas.Trans y = alpha * Aᴴ * x + beta * y if trans = blas.ConjTrans
where alpha and beta are scalars, x and y are vectors, and A is an m×n band matrix with kL sub-diagonals and kU super-diagonals.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cgemm ¶
func (Implementation) Cgemm(tA, tB blas.Transpose, m, n, k int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)
Cgemm 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ᵀ or op(X) = Xᴴ,
alpha and beta are scalars, and A, B and C are matrices, with op(A) an m×k matrix, op(B) a k×n matrix and C an m×n matrix.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cgemv ¶
func (Implementation) Cgemv(trans blas.Transpose, m, n int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)
Cgemv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if trans = blas.NoTrans y = alpha * Aᵀ * x + beta * y if trans = blas.Trans y = alpha * Aᴴ * x + beta * y if trans = blas.ConjTrans
where alpha and beta are scalars, x and y are vectors, and A is an m×n dense matrix.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cgerc ¶
func (Implementation) Cgerc(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)
Cgerc performs the rank-one operation
A += alpha * x * yᴴ
where A is an m×n dense matrix, alpha is a scalar, x is an m element vector, and y is an n element vector.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cgeru ¶
func (Implementation) Cgeru(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)
Cgeru performs the rank-one operation
A += alpha * x * yᵀ
where A is an m×n dense matrix, alpha is a scalar, x is an m element vector, and y is an n element vector.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chbmv ¶
func (Implementation) Chbmv(uplo blas.Uplo, n, k int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)
Chbmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian band matrix with k super-diagonals. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chemm ¶
func (Implementation) Chemm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)
Chemm performs one of the matrix-matrix operations
C = alpha*A*B + beta*C if side == blas.Left C = alpha*B*A + beta*C if side == blas.Right
where alpha and beta are scalars, A is an m×m or n×n hermitian matrix and B and C are m×n matrices. The imaginary parts of the diagonal elements of A are assumed to be zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chemv ¶
func (Implementation) Chemv(uplo blas.Uplo, n int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)
Chemv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian matrix. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cher ¶
func (Implementation) Cher(uplo blas.Uplo, n int, alpha float32, x []complex64, incX int, a []complex64, lda int)
Cher performs the Hermitian rank-one operation
A += alpha * x * xᴴ
where A is an n×n Hermitian matrix, alpha is a real scalar, and x is an n element vector. On entry, the imaginary parts of the diagonal elements of A are ignored and assumed to be zero, on return they will be set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cher2 ¶
func (Implementation) Cher2(uplo blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)
Cher2 performs the Hermitian rank-two operation
A += alpha * x * yᴴ + conj(alpha) * y * xᴴ
where alpha is a scalar, x and y are n element vectors and A is an n×n Hermitian matrix. On entry, the imaginary parts of the diagonal elements are ignored and assumed to be zero. On return they will be set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cher2k ¶
func (Implementation) Cher2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta float32, c []complex64, ldc int)
Cher2k performs one of the hermitian rank-2k operations
C = alpha*A*Bᴴ + conj(alpha)*B*Aᴴ + beta*C if trans == blas.NoTrans C = alpha*Aᴴ*B + conj(alpha)*Bᴴ*A + beta*C if trans == blas.ConjTrans
where alpha and beta are scalars with beta real, C is an n×n hermitian matrix and A and B are n×k matrices in the first case and k×n matrices in the second case.
The imaginary parts of the diagonal elements of C are assumed to be zero, and on return they will be set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cherk ¶
func (Implementation) Cherk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha float32, a []complex64, lda int, beta float32, c []complex64, ldc int)
Cherk performs one of the hermitian rank-k operations
C = alpha*A*Aᴴ + beta*C if trans == blas.NoTrans C = alpha*Aᴴ*A + beta*C if trans == blas.ConjTrans
where alpha and beta are real scalars, C is an n×n hermitian matrix and A is an n×k matrix in the first case and a k×n matrix in the second case.
The imaginary parts of the diagonal elements of C are assumed to be zero, and on return they will be set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chpmv ¶
func (Implementation) Chpmv(uplo blas.Uplo, n int, alpha complex64, ap []complex64, x []complex64, incX int, beta complex64, y []complex64, incY int)
Chpmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian matrix in packed form. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chpr ¶
func (Implementation) Chpr(uplo blas.Uplo, n int, alpha float32, x []complex64, incX int, ap []complex64)
Chpr performs the Hermitian rank-1 operation
A += alpha * x * xᴴ
where alpha is a real scalar, x is a vector, and A is an n×n hermitian matrix in packed form. On entry, the imaginary parts of the diagonal elements are assumed to be zero, and on return they are set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Chpr2 ¶
func (Implementation) Chpr2(uplo blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, ap []complex64)
Chpr2 performs the Hermitian rank-2 operation
A += alpha * x * yᴴ + conj(alpha) * y * xᴴ
where alpha is a complex scalar, x and y are n element vectors, and A is an n×n Hermitian matrix, supplied in packed form. On entry, the imaginary parts of the diagonal elements are assumed to be zero, and on return they are set to zero.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cscal ¶
func (Implementation) Cscal(n int, alpha complex64, x []complex64, incX int)
Cscal scales the vector x by a complex scalar alpha. Cscal has no effect if incX < 0.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Csscal ¶
func (Implementation) Csscal(n int, alpha float32, x []complex64, incX int)
Csscal scales the vector x by a real scalar alpha. Csscal has no effect if incX < 0.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Cswap ¶
Cswap exchanges the elements of two complex vectors x and y.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Csymm ¶
func (Implementation) Csymm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)
Csymm performs one of the matrix-matrix operations
C = alpha*A*B + beta*C if side == blas.Left C = alpha*B*A + beta*C if side == blas.Right
where alpha and beta are scalars, A is an m×m or n×n symmetric matrix and B and C are m×n matrices.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Csyr2k ¶
func (Implementation) Csyr2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)
Csyr2k performs one of the symmetric rank-2k operations
C = alpha*A*Bᵀ + alpha*B*Aᵀ + beta*C if trans == blas.NoTrans C = alpha*Aᵀ*B + alpha*Bᵀ*A + beta*C if trans == blas.Trans
where alpha and beta are scalars, C is an n×n symmetric matrix and A and B are n×k matrices in the first case and k×n matrices in the second case.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Csyrk ¶
func (Implementation) Csyrk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, beta complex64, c []complex64, ldc int)
Csyrk performs one of the symmetric rank-k operations
C = alpha*A*Aᵀ + beta*C if trans == blas.NoTrans C = alpha*Aᵀ*A + beta*C if trans == blas.Trans
where alpha and beta are scalars, C is an n×n symmetric matrix and A is an n×k matrix in the first case and a k×n matrix in the second case.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctbmv ¶
func (Implementation) Ctbmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex64, lda int, x []complex64, incX int)
Ctbmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is an n element vector and A is an n×n triangular band matrix, with (k+1) diagonals.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctbsv ¶
func (Implementation) Ctbsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex64, lda int, x []complex64, incX int)
Ctbsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular band matrix with (k+1) diagonals.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctpmv ¶
func (Implementation) Ctpmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex64, x []complex64, incX int)
Ctpmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is an n element vector and A is an n×n triangular matrix, supplied in packed form.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctpsv ¶
func (Implementation) Ctpsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex64, x []complex64, incX int)
Ctpsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular matrix in packed form.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctrmm ¶
func (Implementation) Ctrmm(side blas.Side, uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int)
Ctrmm performs one of the matrix-matrix operations
B = alpha * op(A) * B if side == blas.Left, B = alpha * B * op(A) if side == blas.Right,
where alpha is a scalar, B is an m×n matrix, A is a unit, or non-unit, upper or lower triangular matrix and op(A) is one of
op(A) = A if trans == blas.NoTrans, op(A) = Aᵀ if trans == blas.Trans, op(A) = Aᴴ if trans == blas.ConjTrans.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctrmv ¶
func (Implementation) Ctrmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex64, lda int, x []complex64, incX int)
Ctrmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is a vector, and A is an n×n triangular matrix.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctrsm ¶
func (Implementation) Ctrsm(side blas.Side, uplo blas.Uplo, transA blas.Transpose, diag blas.Diag, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int)
Ctrsm solves one of the matrix equations
op(A) * X = alpha * B if side == blas.Left, X * op(A) = alpha * B if side == blas.Right,
where alpha is a scalar, X and B are m×n matrices, A is a unit or non-unit, upper or lower triangular matrix and op(A) is one of
op(A) = A if transA == blas.NoTrans, op(A) = Aᵀ if transA == blas.Trans, op(A) = Aᴴ if transA == blas.ConjTrans.
On return the matrix X is overwritten on B.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Ctrsv ¶
func (Implementation) Ctrsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex64, lda int, x []complex64, incX int)
Ctrsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular matrix.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Dasum ¶
func (Implementation) Dasum(n int, x []float64, incX int) float64
Dasum computes the sum of the absolute values of the elements of x.
\sum_i |x[i]|
Dasum returns 0 if incX is negative.
func (Implementation) Dcopy ¶
Dcopy copies the elements of x into the elements of y.
y[i] = x[i] for all i
func (Implementation) Dgbmv ¶
func (Implementation) Dgbmv(tA blas.Transpose, m, n, kL, kU int, alpha float64, a []float64, lda int, x []float64, incX int, beta float64, y []float64, incY int)
Dgbmv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if tA == blas.NoTrans y = alpha * Aᵀ * x + beta * y if tA == blas.Trans or blas.ConjTrans
where A is an m×n band matrix with kL sub-diagonals and kU super-diagonals, x and y are vectors, and alpha and beta are scalars.
func (Implementation) Dgemm ¶
func (Implementation) Dgemm(tA, tB blas.Transpose, m, n, k int, alpha float64, a []float64, lda int, b []float64, ldb int, beta float64, c []float64, ldc int)
Dgemm performs one of the matrix-matrix operations
C = alpha * A * B + beta * C C = alpha * Aᵀ * B + beta * C C = alpha * A * Bᵀ + beta * C C = alpha * Aᵀ * Bᵀ + beta * C
where A is an m×k or k×m dense matrix, B is an n×k or k×n dense matrix, C is an m×n matrix, and alpha and beta are scalars. tA and tB specify whether A or B are transposed.
func (Implementation) Dgemv ¶
func (Implementation) Dgemv(tA blas.Transpose, m, n int, alpha float64, a []float64, lda int, x []float64, incX int, beta float64, y []float64, incY int)
Dgemv computes
y = alpha * A * x + beta * y if tA = blas.NoTrans y = alpha * Aᵀ * x + beta * y if tA = blas.Trans or blas.ConjTrans
where A is an m×n dense matrix, x and y are vectors, and alpha and beta are scalars.
func (Implementation) Dger ¶
func (Implementation) Dger(m, n int, alpha float64, x []float64, incX int, y []float64, incY int, a []float64, lda int)
Dger performs the rank-one operation
A += alpha * x * yᵀ
where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.
func (Implementation) Dnrm2 ¶
func (Implementation) Dnrm2(n int, x []float64, incX int) float64
Dnrm2 computes the Euclidean norm of a vector,
sqrt(\sum_i x[i] * x[i]).
This function returns 0 if incX is negative.
func (Implementation) Drot ¶
func (Implementation) Drot(n int, x []float64, incX int, y []float64, incY int, c float64, s float64)
Drot applies a plane transformation.
x[i] = c * x[i] + s * y[i] y[i] = c * y[i] - s * x[i]
func (Implementation) Drotg ¶
func (Implementation) Drotg(a, b float64) (c, s, r, z float64)
Drotg computes the plane rotation
_ _ _ _ _ _ | c s | | a | | r | | -s c | * | b | = | 0 | ‾ ‾ ‾ ‾ ‾ ‾
where
r = ±√(a^2 + b^2) c = a/r, the cosine of the plane rotation s = b/r, the sine of the plane rotation
NOTE: There is a discrepancy between the reference implementation and the BLAS technical manual regarding the sign for r when a or b are zero. Drotg agrees with the definition in the manual and other common BLAS implementations.
func (Implementation) Drotm ¶
func (Implementation) Drotm(n int, x []float64, incX int, y []float64, incY int, p blas.DrotmParams)
Drotm applies the modified Givens rotation to the 2×n matrix.
func (Implementation) Drotmg ¶
func (Implementation) Drotmg(d1, d2, x1, y1 float64) (p blas.DrotmParams, rd1, rd2, rx1 float64)
Drotmg computes the modified Givens rotation. See http://www.netlib.org/lapack/explore-html/df/deb/drotmg_8f.html for more details.
func (Implementation) Dsbmv ¶
func (Implementation) Dsbmv(ul blas.Uplo, n, k int, alpha float64, a []float64, lda int, x []float64, incX int, beta float64, y []float64, incY int)
Dsbmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric band matrix with k super-diagonals, x and y are vectors, and alpha and beta are scalars.
func (Implementation) Dscal ¶
func (Implementation) Dscal(n int, alpha float64, x []float64, incX int)
Dscal scales x by alpha.
x[i] *= alpha
Dscal has no effect if incX < 0.
func (Implementation) Dsdot ¶
Dsdot computes the dot product of the two vectors
\sum_i x[i]*y[i]
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Dspmv ¶
func (Implementation) Dspmv(ul blas.Uplo, n int, alpha float64, ap []float64, x []float64, incX int, beta float64, y []float64, incY int)
Dspmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha and beta are scalars.
func (Implementation) Dspr ¶
Dspr performs the symmetric rank-one operation
A += alpha * x * xᵀ
where A is an n×n symmetric matrix in packed format, x is a vector, and alpha is a scalar.
func (Implementation) Dspr2 ¶
func (Implementation) Dspr2(ul blas.Uplo, n int, alpha float64, x []float64, incX int, y []float64, incY int, ap []float64)
Dspr2 performs the symmetric rank-2 update
A += alpha * x * yᵀ + alpha * y * xᵀ
where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha is a scalar.
func (Implementation) Dswap ¶
Dswap exchanges the elements of two vectors.
x[i], y[i] = y[i], x[i] for all i
func (Implementation) Dsymm ¶
func (Implementation) Dsymm(s blas.Side, ul blas.Uplo, m, n int, alpha float64, a []float64, lda int, b []float64, ldb int, beta float64, c []float64, ldc int)
Dsymm performs one of the matrix-matrix operations
C = alpha * A * B + beta * C if side == blas.Left C = alpha * B * A + beta * C if side == blas.Right
where A is an n×n or m×m symmetric matrix, B and C are m×n matrices, and alpha is a scalar.
func (Implementation) Dsymv ¶
func (Implementation) Dsymv(ul blas.Uplo, n int, alpha float64, a []float64, lda int, x []float64, incX int, beta float64, y []float64, incY int)
Dsymv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric matrix, x and y are vectors, and alpha and beta are scalars.
func (Implementation) Dsyr ¶
func (Implementation) Dsyr(ul blas.Uplo, n int, alpha float64, x []float64, incX int, a []float64, lda int)
Dsyr performs the symmetric rank-one update
A += alpha * x * xᵀ
where A is an n×n symmetric matrix, and x is a vector.
func (Implementation) Dsyr2 ¶
func (Implementation) Dsyr2(ul blas.Uplo, n int, alpha float64, x []float64, incX int, y []float64, incY int, a []float64, lda int)
Dsyr2 performs the symmetric rank-two update
A += alpha * x * yᵀ + alpha * y * xᵀ
where A is an n×n symmetric matrix, x and y are vectors, and alpha is a scalar.
func (Implementation) Dsyr2k ¶
func (Implementation) Dsyr2k(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float64, a []float64, lda int, b []float64, ldb int, beta float64, c []float64, ldc int)
Dsyr2k performs one of the symmetric rank 2k operations
C = alpha * A * Bᵀ + alpha * B * Aᵀ + beta * C if tA == blas.NoTrans C = alpha * Aᵀ * B + alpha * Bᵀ * A + beta * C if tA == blas.Trans or tA == blas.ConjTrans
where A and B are n×k or k×n matrices, C is an n×n symmetric matrix, and alpha and beta are scalars.
func (Implementation) Dsyrk ¶
func (Implementation) Dsyrk(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float64, a []float64, lda int, beta float64, c []float64, ldc int)
Dsyrk performs one of the symmetric rank-k operations
C = alpha * A * Aᵀ + beta * C if tA == blas.NoTrans C = alpha * Aᵀ * A + beta * C if tA == blas.Trans or tA == blas.ConjTrans
where A is an n×k or k×n matrix, C is an n×n symmetric matrix, and alpha and beta are scalars.
func (Implementation) Dtbmv ¶
func (Implementation) Dtbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float64, lda int, x []float64, incX int)
Dtbmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular band matrix with k+1 diagonals, and x is a vector.
func (Implementation) Dtbsv ¶
func (Implementation) Dtbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float64, lda int, x []float64, incX int)
Dtbsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or tA == blas.ConjTrans
where A is an n×n triangular band matrix with k+1 diagonals, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
func (Implementation) Dtpmv ¶
func (Implementation) Dtpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float64, x []float64, incX int)
Dtpmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix in packed format, and x is a vector.
func (Implementation) Dtpsv ¶
func (Implementation) Dtpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float64, x []float64, incX int)
Dtpsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix in packed format, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
func (Implementation) Dtrmm ¶
func (Implementation) Dtrmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha float64, a []float64, lda int, b []float64, ldb int)
Dtrmm performs one of the matrix-matrix operations
B = alpha * A * B if tA == blas.NoTrans and side == blas.Left B = alpha * Aᵀ * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Left B = alpha * B * A if tA == blas.NoTrans and side == blas.Right B = alpha * B * Aᵀ if tA == blas.Trans or blas.ConjTrans, and side == blas.Right
where A is an n×n or m×m triangular matrix, B is an m×n matrix, and alpha is a scalar.
func (Implementation) Dtrmv ¶
func (Implementation) Dtrmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float64, lda int, x []float64, incX int)
Dtrmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix, and x is a vector.
func (Implementation) Dtrsm ¶
func (Implementation) Dtrsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha float64, a []float64, lda int, b []float64, ldb int)
Dtrsm solves one of the matrix equations
A * X = alpha * B if tA == blas.NoTrans and side == blas.Left Aᵀ * X = alpha * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Left X * A = alpha * B if tA == blas.NoTrans and side == blas.Right X * Aᵀ = alpha * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Right
where A is an n×n or m×m triangular matrix, X and B are m×n matrices, and alpha is a scalar.
At entry to the function, X contains the values of B, and the result is stored in-place into X.
No check is made that A is invertible.
func (Implementation) Dtrsv ¶
func (Implementation) Dtrsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float64, lda int, x []float64, incX int)
Dtrsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
func (Implementation) Dzasum ¶
func (Implementation) Dzasum(n int, x []complex128, incX int) float64
Dzasum returns the sum of the absolute values of the elements of x
\sum_i |Re(x[i])| + |Im(x[i])|
Dzasum returns 0 if incX is negative.
func (Implementation) Dznrm2 ¶
func (Implementation) Dznrm2(n int, x []complex128, incX int) float64
Dznrm2 computes the Euclidean norm of the complex vector x,
‖x‖_2 = sqrt(\sum_i x[i] * conj(x[i])).
This function returns 0 if incX is negative.
func (Implementation) Icamax ¶
func (Implementation) Icamax(n int, x []complex64, incX int) int
Icamax returns the index of the first element of x having largest |Re(·)|+|Im(·)|. Icamax returns -1 if n is 0 or incX is negative.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Idamax ¶
func (Implementation) Idamax(n int, x []float64, incX int) int
Idamax returns the index of an element of x with the largest absolute value. If there are multiple such indices the earliest is returned. Idamax returns -1 if n == 0.
func (Implementation) Isamax ¶
func (Implementation) Isamax(n int, x []float32, incX int) int
Isamax returns the index of an element of x with the largest absolute value. If there are multiple such indices the earliest is returned. Isamax returns -1 if n == 0.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Izamax ¶
func (Implementation) Izamax(n int, x []complex128, incX int) int
Izamax returns the index of the first element of x having largest |Re(·)|+|Im(·)|. Izamax returns -1 if n is 0 or incX is negative.
func (Implementation) Sasum ¶
func (Implementation) Sasum(n int, x []float32, incX int) float32
Sasum computes the sum of the absolute values of the elements of x.
\sum_i |x[i]|
Sasum returns 0 if incX is negative.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Saxpy ¶
Saxpy adds alpha times x to y
y[i] += alpha * x[i] for all i
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Scasum ¶
func (Implementation) Scasum(n int, x []complex64, incX int) float32
Scasum returns the sum of the absolute values of the elements of x
\sum_i |Re(x[i])| + |Im(x[i])|
Scasum returns 0 if incX is negative.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Scnrm2 ¶
func (Implementation) Scnrm2(n int, x []complex64, incX int) float32
Scnrm2 computes the Euclidean norm of the complex vector x,
‖x‖_2 = sqrt(\sum_i x[i] * conj(x[i])).
This function returns 0 if incX is negative.
Complex64 implementations are autogenerated and not directly tested.
func (Implementation) Scopy ¶
Scopy copies the elements of x into the elements of y.
y[i] = x[i] for all i
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sdot ¶
Sdot computes the dot product of the two vectors
\sum_i x[i]*y[i]
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sdsdot ¶
func (Implementation) Sdsdot(n int, alpha float32, x []float32, incX int, y []float32, incY int) float32
Sdsdot computes the dot product of the two vectors plus a constant
alpha + \sum_i x[i]*y[i]
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sgbmv ¶
func (Implementation) Sgbmv(tA blas.Transpose, m, n, kL, kU int, alpha float32, a []float32, lda int, x []float32, incX int, beta float32, y []float32, incY int)
Sgbmv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if tA == blas.NoTrans y = alpha * Aᵀ * x + beta * y if tA == blas.Trans or blas.ConjTrans
where A is an m×n band matrix with kL sub-diagonals and kU super-diagonals, x and y are vectors, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sgemm ¶
func (Implementation) Sgemm(tA, tB blas.Transpose, m, n, k int, alpha float32, a []float32, lda int, b []float32, ldb int, beta float32, c []float32, ldc int)
Sgemm performs one of the matrix-matrix operations
C = alpha * A * B + beta * C C = alpha * Aᵀ * B + beta * C C = alpha * A * Bᵀ + beta * C C = alpha * Aᵀ * Bᵀ + beta * C
where A is an m×k or k×m dense matrix, B is an n×k or k×n dense matrix, C is an m×n matrix, and alpha and beta are scalars. tA and tB specify whether A or B are transposed.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sgemv ¶
func (Implementation) Sgemv(tA blas.Transpose, m, n int, alpha float32, a []float32, lda int, x []float32, incX int, beta float32, y []float32, incY int)
Sgemv computes
y = alpha * A * x + beta * y if tA = blas.NoTrans y = alpha * Aᵀ * x + beta * y if tA = blas.Trans or blas.ConjTrans
where A is an m×n dense matrix, x and y are vectors, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sger ¶
func (Implementation) Sger(m, n int, alpha float32, x []float32, incX int, y []float32, incY int, a []float32, lda int)
Sger performs the rank-one operation
A += alpha * x * yᵀ
where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Snrm2 ¶
func (Implementation) Snrm2(n int, x []float32, incX int) float32
Snrm2 computes the Euclidean norm of a vector,
sqrt(\sum_i x[i] * x[i]).
This function returns 0 if incX is negative.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Srot ¶
func (Implementation) Srot(n int, x []float32, incX int, y []float32, incY int, c float32, s float32)
Srot applies a plane transformation.
x[i] = c * x[i] + s * y[i] y[i] = c * y[i] - s * x[i]
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Srotg ¶
func (Implementation) Srotg(a, b float32) (c, s, r, z float32)
Srotg computes the plane rotation
_ _ _ _ _ _ | c s | | a | | r | | -s c | * | b | = | 0 | ‾ ‾ ‾ ‾ ‾ ‾
where
r = ±√(a^2 + b^2) c = a/r, the cosine of the plane rotation s = b/r, the sine of the plane rotation
NOTE: There is a discrepancy between the reference implementation and the BLAS technical manual regarding the sign for r when a or b are zero. Srotg agrees with the definition in the manual and other common BLAS implementations.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Srotm ¶
func (Implementation) Srotm(n int, x []float32, incX int, y []float32, incY int, p blas.SrotmParams)
Srotm applies the modified Givens rotation to the 2×n matrix.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Srotmg ¶
func (Implementation) Srotmg(d1, d2, x1, y1 float32) (p blas.SrotmParams, rd1, rd2, rx1 float32)
Srotmg computes the modified Givens rotation. See http://www.netlib.org/lapack/explore-html/df/deb/drotmg_8f.html for more details.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssbmv ¶
func (Implementation) Ssbmv(ul blas.Uplo, n, k int, alpha float32, a []float32, lda int, x []float32, incX int, beta float32, y []float32, incY int)
Ssbmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric band matrix with k super-diagonals, x and y are vectors, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sscal ¶
func (Implementation) Sscal(n int, alpha float32, x []float32, incX int)
Sscal scales x by alpha.
x[i] *= alpha
Sscal has no effect if incX < 0.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sspmv ¶
func (Implementation) Sspmv(ul blas.Uplo, n int, alpha float32, ap []float32, x []float32, incX int, beta float32, y []float32, incY int)
Sspmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sspr ¶
Sspr performs the symmetric rank-one operation
A += alpha * x * xᵀ
where A is an n×n symmetric matrix in packed format, x is a vector, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sspr2 ¶
func (Implementation) Sspr2(ul blas.Uplo, n int, alpha float32, x []float32, incX int, y []float32, incY int, ap []float32)
Sspr2 performs the symmetric rank-2 update
A += alpha * x * yᵀ + alpha * y * xᵀ
where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Sswap ¶
Sswap exchanges the elements of two vectors.
x[i], y[i] = y[i], x[i] for all i
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssymm ¶
func (Implementation) Ssymm(s blas.Side, ul blas.Uplo, m, n int, alpha float32, a []float32, lda int, b []float32, ldb int, beta float32, c []float32, ldc int)
Ssymm performs one of the matrix-matrix operations
C = alpha * A * B + beta * C if side == blas.Left C = alpha * B * A + beta * C if side == blas.Right
where A is an n×n or m×m symmetric matrix, B and C are m×n matrices, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssymv ¶
func (Implementation) Ssymv(ul blas.Uplo, n int, alpha float32, a []float32, lda int, x []float32, incX int, beta float32, y []float32, incY int)
Ssymv performs the matrix-vector operation
y = alpha * A * x + beta * y
where A is an n×n symmetric matrix, x and y are vectors, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssyr ¶
func (Implementation) Ssyr(ul blas.Uplo, n int, alpha float32, x []float32, incX int, a []float32, lda int)
Ssyr performs the symmetric rank-one update
A += alpha * x * xᵀ
where A is an n×n symmetric matrix, and x is a vector.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssyr2 ¶
func (Implementation) Ssyr2(ul blas.Uplo, n int, alpha float32, x []float32, incX int, y []float32, incY int, a []float32, lda int)
Ssyr2 performs the symmetric rank-two update
A += alpha * x * yᵀ + alpha * y * xᵀ
where A is an n×n symmetric matrix, x and y are vectors, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssyr2k ¶
func (Implementation) Ssyr2k(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float32, a []float32, lda int, b []float32, ldb int, beta float32, c []float32, ldc int)
Ssyr2k performs one of the symmetric rank 2k operations
C = alpha * A * Bᵀ + alpha * B * Aᵀ + beta * C if tA == blas.NoTrans C = alpha * Aᵀ * B + alpha * Bᵀ * A + beta * C if tA == blas.Trans or tA == blas.ConjTrans
where A and B are n×k or k×n matrices, C is an n×n symmetric matrix, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Ssyrk ¶
func (Implementation) Ssyrk(ul blas.Uplo, tA blas.Transpose, n, k int, alpha float32, a []float32, lda int, beta float32, c []float32, ldc int)
Ssyrk performs one of the symmetric rank-k operations
C = alpha * A * Aᵀ + beta * C if tA == blas.NoTrans C = alpha * Aᵀ * A + beta * C if tA == blas.Trans or tA == blas.ConjTrans
where A is an n×k or k×n matrix, C is an n×n symmetric matrix, and alpha and beta are scalars.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Stbmv ¶
func (Implementation) Stbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float32, lda int, x []float32, incX int)
Stbmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular band matrix with k+1 diagonals, and x is a vector.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Stbsv ¶
func (Implementation) Stbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []float32, lda int, x []float32, incX int)
Stbsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or tA == blas.ConjTrans
where A is an n×n triangular band matrix with k+1 diagonals, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Stpmv ¶
func (Implementation) Stpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float32, x []float32, incX int)
Stpmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix in packed format, and x is a vector.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Stpsv ¶
func (Implementation) Stpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap []float32, x []float32, incX int)
Stpsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix in packed format, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Strmm ¶
func (Implementation) Strmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha float32, a []float32, lda int, b []float32, ldb int)
Strmm performs one of the matrix-matrix operations
B = alpha * A * B if tA == blas.NoTrans and side == blas.Left B = alpha * Aᵀ * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Left B = alpha * B * A if tA == blas.NoTrans and side == blas.Right B = alpha * B * Aᵀ if tA == blas.Trans or blas.ConjTrans, and side == blas.Right
where A is an n×n or m×m triangular matrix, B is an m×n matrix, and alpha is a scalar.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Strmv ¶
func (Implementation) Strmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float32, lda int, x []float32, incX int)
Strmv performs one of the matrix-vector operations
x = A * x if tA == blas.NoTrans x = Aᵀ * x if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix, and x is a vector.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Strsm ¶
func (Implementation) Strsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha float32, a []float32, lda int, b []float32, ldb int)
Strsm solves one of the matrix equations
A * X = alpha * B if tA == blas.NoTrans and side == blas.Left Aᵀ * X = alpha * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Left X * A = alpha * B if tA == blas.NoTrans and side == blas.Right X * Aᵀ = alpha * B if tA == blas.Trans or blas.ConjTrans, and side == blas.Right
where A is an n×n or m×m triangular matrix, X and B are m×n matrices, and alpha is a scalar.
At entry to the function, X contains the values of B, and the result is stored in-place into X.
No check is made that A is invertible.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Strsv ¶
func (Implementation) Strsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []float32, lda int, x []float32, incX int)
Strsv solves one of the systems of equations
A * x = b if tA == blas.NoTrans Aᵀ * x = b if tA == blas.Trans or blas.ConjTrans
where A is an n×n triangular matrix, and x and b are vectors.
At entry to the function, x contains the values of b, and the result is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Float32 implementations are autogenerated and not directly tested.
func (Implementation) Zaxpy ¶
func (Implementation) Zaxpy(n int, alpha complex128, x []complex128, incX int, y []complex128, incY int)
Zaxpy adds alpha times x to y:
y[i] += alpha * x[i] for all i
func (Implementation) Zcopy ¶
func (Implementation) Zcopy(n int, x []complex128, incX int, y []complex128, incY int)
Zcopy copies the vector x to vector y.
func (Implementation) Zdotc ¶
func (Implementation) Zdotc(n int, x []complex128, incX int, y []complex128, incY int) complex128
Zdotc computes the dot product
xᴴ · y
of two complex vectors x and y.
func (Implementation) Zdotu ¶
func (Implementation) Zdotu(n int, x []complex128, incX int, y []complex128, incY int) complex128
Zdotu computes the dot product
xᵀ · y
of two complex vectors x and y.
func (Implementation) Zdscal ¶
func (Implementation) Zdscal(n int, alpha float64, x []complex128, incX int)
Zdscal scales the vector x by a real scalar alpha. Zdscal has no effect if incX < 0.
func (Implementation) Zgbmv ¶
func (Implementation) Zgbmv(trans blas.Transpose, m, n, kL, kU int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)
Zgbmv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if trans = blas.NoTrans y = alpha * Aᵀ * x + beta * y if trans = blas.Trans y = alpha * Aᴴ * x + beta * y if trans = blas.ConjTrans
where alpha and beta are scalars, x and y are vectors, and A is an m×n band matrix with kL sub-diagonals and kU super-diagonals.
func (Implementation) Zgemm ¶
func (Implementation) Zgemm(tA, tB blas.Transpose, m, n, k int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)
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ᵀ or op(X) = Xᴴ,
alpha and beta are scalars, and A, B and C are matrices, with op(A) an m×k matrix, op(B) a k×n matrix and C an m×n matrix.
func (Implementation) Zgemv ¶
func (Implementation) Zgemv(trans blas.Transpose, m, n int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)
Zgemv performs one of the matrix-vector operations
y = alpha * A * x + beta * y if trans = blas.NoTrans y = alpha * Aᵀ * x + beta * y if trans = blas.Trans y = alpha * Aᴴ * x + beta * y if trans = blas.ConjTrans
where alpha and beta are scalars, x and y are vectors, and A is an m×n dense matrix.
func (Implementation) Zgerc ¶
func (Implementation) Zgerc(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)
Zgerc performs the rank-one operation
A += alpha * x * yᴴ
where A is an m×n dense matrix, alpha is a scalar, x is an m element vector, and y is an n element vector.
func (Implementation) Zgeru ¶
func (Implementation) Zgeru(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)
Zgeru performs the rank-one operation
A += alpha * x * yᵀ
where A is an m×n dense matrix, alpha is a scalar, x is an m element vector, and y is an n element vector.
func (Implementation) Zhbmv ¶
func (Implementation) Zhbmv(uplo blas.Uplo, n, k int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)
Zhbmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian band matrix with k super-diagonals. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
func (Implementation) Zhemm ¶
func (Implementation) Zhemm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)
Zhemm performs one of the matrix-matrix operations
C = alpha*A*B + beta*C if side == blas.Left C = alpha*B*A + beta*C if side == blas.Right
where alpha and beta are scalars, A is an m×m or n×n hermitian matrix and B and C are m×n matrices. The imaginary parts of the diagonal elements of A are assumed to be zero.
func (Implementation) Zhemv ¶
func (Implementation) Zhemv(uplo blas.Uplo, n int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)
Zhemv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian matrix. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
func (Implementation) Zher ¶
func (Implementation) Zher(uplo blas.Uplo, n int, alpha float64, x []complex128, incX int, a []complex128, lda int)
Zher performs the Hermitian rank-one operation
A += alpha * x * xᴴ
where A is an n×n Hermitian matrix, alpha is a real scalar, and x is an n element vector. On entry, the imaginary parts of the diagonal elements of A are ignored and assumed to be zero, on return they will be set to zero.
func (Implementation) Zher2 ¶
func (Implementation) Zher2(uplo blas.Uplo, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)
Zher2 performs the Hermitian rank-two operation
A += alpha * x * yᴴ + conj(alpha) * y * xᴴ
where alpha is a scalar, x and y are n element vectors and A is an n×n Hermitian matrix. On entry, the imaginary parts of the diagonal elements are ignored and assumed to be zero. On return they will be set to zero.
func (Implementation) Zher2k ¶
func (Implementation) Zher2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta float64, c []complex128, ldc int)
Zher2k performs one of the hermitian rank-2k operations
C = alpha*A*Bᴴ + conj(alpha)*B*Aᴴ + beta*C if trans == blas.NoTrans C = alpha*Aᴴ*B + conj(alpha)*Bᴴ*A + beta*C if trans == blas.ConjTrans
where alpha and beta are scalars with beta real, C is an n×n hermitian matrix and A and B are n×k matrices in the first case and k×n matrices in the second case.
The imaginary parts of the diagonal elements of C are assumed to be zero, and on return they will be set to zero.
func (Implementation) Zherk ¶
func (Implementation) Zherk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha float64, a []complex128, lda int, beta float64, c []complex128, ldc int)
Zherk performs one of the hermitian rank-k operations
C = alpha*A*Aᴴ + beta*C if trans == blas.NoTrans C = alpha*Aᴴ*A + beta*C if trans == blas.ConjTrans
where alpha and beta are real scalars, C is an n×n hermitian matrix and A is an n×k matrix in the first case and a k×n matrix in the second case.
The imaginary parts of the diagonal elements of C are assumed to be zero, and on return they will be set to zero.
func (Implementation) Zhpmv ¶
func (Implementation) Zhpmv(uplo blas.Uplo, n int, alpha complex128, ap []complex128, x []complex128, incX int, beta complex128, y []complex128, incY int)
Zhpmv performs the matrix-vector operation
y = alpha * A * x + beta * y
where alpha and beta are scalars, x and y are vectors, and A is an n×n Hermitian matrix in packed form. The imaginary parts of the diagonal elements of A are ignored and assumed to be zero.
func (Implementation) Zhpr ¶
func (Implementation) Zhpr(uplo blas.Uplo, n int, alpha float64, x []complex128, incX int, ap []complex128)
Zhpr performs the Hermitian rank-1 operation
A += alpha * x * xᴴ
where alpha is a real scalar, x is a vector, and A is an n×n hermitian matrix in packed form. On entry, the imaginary parts of the diagonal elements are assumed to be zero, and on return they are set to zero.
func (Implementation) Zhpr2 ¶
func (Implementation) Zhpr2(uplo blas.Uplo, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, ap []complex128)
Zhpr2 performs the Hermitian rank-2 operation
A += alpha * x * yᴴ + conj(alpha) * y * xᴴ
where alpha is a complex scalar, x and y are n element vectors, and A is an n×n Hermitian matrix, supplied in packed form. On entry, the imaginary parts of the diagonal elements are assumed to be zero, and on return they are set to zero.
func (Implementation) Zscal ¶
func (Implementation) Zscal(n int, alpha complex128, x []complex128, incX int)
Zscal scales the vector x by a complex scalar alpha. Zscal has no effect if incX < 0.
func (Implementation) Zswap ¶
func (Implementation) Zswap(n int, x []complex128, incX int, y []complex128, incY int)
Zswap exchanges the elements of two complex vectors x and y.
func (Implementation) Zsymm ¶
func (Implementation) Zsymm(side blas.Side, uplo blas.Uplo, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)
Zsymm performs one of the matrix-matrix operations
C = alpha*A*B + beta*C if side == blas.Left C = alpha*B*A + beta*C if side == blas.Right
where alpha and beta are scalars, A is an m×m or n×n symmetric matrix and B and C are m×n matrices.
func (Implementation) Zsyr2k ¶
func (Implementation) Zsyr2k(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)
Zsyr2k performs one of the symmetric rank-2k operations
C = alpha*A*Bᵀ + alpha*B*Aᵀ + beta*C if trans == blas.NoTrans C = alpha*Aᵀ*B + alpha*Bᵀ*A + beta*C if trans == blas.Trans
where alpha and beta are scalars, C is an n×n symmetric matrix and A and B are n×k matrices in the first case and k×n matrices in the second case.
func (Implementation) Zsyrk ¶
func (Implementation) Zsyrk(uplo blas.Uplo, trans blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, beta complex128, c []complex128, ldc int)
Zsyrk performs one of the symmetric rank-k operations
C = alpha*A*Aᵀ + beta*C if trans == blas.NoTrans C = alpha*Aᵀ*A + beta*C if trans == blas.Trans
where alpha and beta are scalars, C is an n×n symmetric matrix and A is an n×k matrix in the first case and a k×n matrix in the second case.
func (Implementation) Ztbmv ¶
func (Implementation) Ztbmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex128, lda int, x []complex128, incX int)
Ztbmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is an n element vector and A is an n×n triangular band matrix, with (k+1) diagonals.
func (Implementation) Ztbsv ¶
func (Implementation) Ztbsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, k int, a []complex128, lda int, x []complex128, incX int)
Ztbsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular band matrix with (k+1) diagonals.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
func (Implementation) Ztpmv ¶
func (Implementation) Ztpmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex128, x []complex128, incX int)
Ztpmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is an n element vector and A is an n×n triangular matrix, supplied in packed form.
func (Implementation) Ztpsv ¶
func (Implementation) Ztpsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, ap []complex128, x []complex128, incX int)
Ztpsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular matrix in packed form.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
func (Implementation) Ztrmm ¶
func (Implementation) Ztrmm(side blas.Side, uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int)
Ztrmm performs one of the matrix-matrix operations
B = alpha * op(A) * B if side == blas.Left, B = alpha * B * op(A) if side == blas.Right,
where alpha is a scalar, B is an m×n matrix, A is a unit, or non-unit, upper or lower triangular matrix and op(A) is one of
op(A) = A if trans == blas.NoTrans, op(A) = Aᵀ if trans == blas.Trans, op(A) = Aᴴ if trans == blas.ConjTrans.
func (Implementation) Ztrmv ¶
func (Implementation) Ztrmv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex128, lda int, x []complex128, incX int)
Ztrmv performs one of the matrix-vector operations
x = A * x if trans = blas.NoTrans x = Aᵀ * x if trans = blas.Trans x = Aᴴ * x if trans = blas.ConjTrans
where x is a vector, and A is an n×n triangular matrix.
func (Implementation) Ztrsm ¶
func (Implementation) Ztrsm(side blas.Side, uplo blas.Uplo, transA blas.Transpose, diag blas.Diag, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int)
Ztrsm solves one of the matrix equations
op(A) * X = alpha * B if side == blas.Left, X * op(A) = alpha * B if side == blas.Right,
where alpha is a scalar, X and B are m×n matrices, A is a unit or non-unit, upper or lower triangular matrix and op(A) is one of
op(A) = A if transA == blas.NoTrans, op(A) = Aᵀ if transA == blas.Trans, op(A) = Aᴴ if transA == blas.ConjTrans.
On return the matrix X is overwritten on B.
func (Implementation) Ztrsv ¶
func (Implementation) Ztrsv(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n int, a []complex128, lda int, x []complex128, incX int)
Ztrsv solves one of the systems of equations
A * x = b if trans == blas.NoTrans Aᵀ * x = b if trans == blas.Trans Aᴴ * x = b if trans == blas.ConjTrans
where b and x are n element vectors and A is an n×n triangular matrix.
On entry, x contains the values of b, and the solution is stored in-place into x.
No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.
Source Files ¶
- dgemm.go
- doc.go
- errors.go
- gemv.go
- gonum.go
- level1cmplx128.go
- level1cmplx64.go
- level1float32.go
- level1float32_dsdot.go
- level1float32_sdot.go
- level1float32_sdsdot.go
- level1float64.go
- level1float64_ddot.go
- level2cmplx128.go
- level2cmplx64.go
- level2float32.go
- level2float64.go
- level3cmplx128.go
- level3cmplx64.go
- level3float32.go
- level3float64.go
- sgemm.go