cgo

package
v0.0.0-...-dfba71b Latest Latest
Warning

This package is not in the latest version of its module.

Go to latest
Published: Dec 9, 2020 License: BSD-2-Clause Imports: 3 Imported by: 0

Documentation

Overview

Package cgo provides bindings to a C BLAS library. This wrapper interface panics when the input arguments are invalid as per the standard, for example if a vector increment is zero. Please note that the treatment of NaN values is not specified, and differs among the BLAS implementations. github.com/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.

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

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

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)

func (Implementation) Ccopy

func (Implementation) Ccopy(n int, x []complex64, incX int, y []complex64, incY int)

func (Implementation) Cdotc

func (Implementation) Cdotc(n int, x []complex64, incX int, y []complex64, incY int) (dotc complex64)

func (Implementation) Cdotu

func (Implementation) Cdotu(n int, x []complex64, incX int, y []complex64, incY int) (dotu complex64)

func (Implementation) Cgbmv

func (Implementation) Cgbmv(tA blas.Transpose, m, n, kL, kU int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)

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)

func (Implementation) Cgemv

func (Implementation) Cgemv(tA blas.Transpose, m, n int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)

func (Implementation) Cgerc

func (Implementation) Cgerc(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)

func (Implementation) Cgeru

func (Implementation) Cgeru(m, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)

func (Implementation) Chbmv

func (Implementation) Chbmv(ul blas.Uplo, n, k int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)

func (Implementation) Chemm

func (Implementation) Chemm(s blas.Side, ul blas.Uplo, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)

func (Implementation) Chemv

func (Implementation) Chemv(ul blas.Uplo, n int, alpha complex64, a []complex64, lda int, x []complex64, incX int, beta complex64, y []complex64, incY int)

func (Implementation) Cher

func (Implementation) Cher(ul blas.Uplo, n int, alpha float32, x []complex64, incX int, a []complex64, lda int)

func (Implementation) Cher2

func (Implementation) Cher2(ul blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, a []complex64, lda int)

func (Implementation) Cher2k

func (Implementation) Cher2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta float32, c []complex64, ldc int)

func (Implementation) Cherk

func (Implementation) Cherk(ul blas.Uplo, t blas.Transpose, n, k int, alpha float32, a []complex64, lda int, beta float32, c []complex64, ldc int)

func (Implementation) Chpmv

func (Implementation) Chpmv(ul blas.Uplo, n int, alpha complex64, ap, x []complex64, incX int, beta complex64, y []complex64, incY int)

func (Implementation) Chpr

func (Implementation) Chpr(ul blas.Uplo, n int, alpha float32, x []complex64, incX int, ap []complex64)

func (Implementation) Chpr2

func (Implementation) Chpr2(ul blas.Uplo, n int, alpha complex64, x []complex64, incX int, y []complex64, incY int, ap []complex64)

func (Implementation) Cscal

func (Implementation) Cscal(n int, alpha complex64, x []complex64, incX int)

func (Implementation) Csscal

func (Implementation) Csscal(n int, alpha float32, x []complex64, incX int)

func (Implementation) Cswap

func (Implementation) Cswap(n int, x []complex64, incX int, y []complex64, incY int)

func (Implementation) Csymm

func (Implementation) Csymm(s blas.Side, ul blas.Uplo, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)

func (Implementation) Csyr2k

func (Implementation) Csyr2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, b []complex64, ldb int, beta complex64, c []complex64, ldc int)

func (Implementation) Csyrk

func (Implementation) Csyrk(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex64, a []complex64, lda int, beta complex64, c []complex64, ldc int)

func (Implementation) Ctbmv

func (Implementation) Ctbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []complex64, lda int, x []complex64, incX int)

func (Implementation) Ctbsv

func (Implementation) Ctbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []complex64, lda int, x []complex64, incX int)

func (Implementation) Ctpmv

func (Implementation) Ctpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []complex64, incX int)

func (Implementation) Ctpsv

func (Implementation) Ctpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []complex64, incX int)

func (Implementation) Ctrmm

func (Implementation) Ctrmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int)

func (Implementation) Ctrmv

func (Implementation) Ctrmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []complex64, lda int, x []complex64, incX int)

func (Implementation) Ctrsm

func (Implementation) Ctrsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha complex64, a []complex64, lda int, b []complex64, ldb int)

func (Implementation) Ctrsv

func (Implementation) Ctrsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []complex64, lda int, x []complex64, incX int)

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

func (Implementation) Daxpy(n int, alpha float64, x []float64, incX int, y []float64, incY int)

Daxpy adds alpha times x to y

y[i] += alpha * x[i] for all i

func (Implementation) Dcopy

func (Implementation) Dcopy(n int, x []float64, incX int, y []float64, incY int)

Dcopy copies the elements of x into the elements of y.

y[i] = x[i] for all i

func (Implementation) Ddot

func (Implementation) Ddot(n int, x []float64, incX int, y []float64, incY int) float64

Ddot computes the dot product of the two vectors

\sum_i x[i]*y[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 computes

y = alpha * A * x + beta * y if tA == blas.NoTrans
y = alpha * A^T * x + beta * y if tA == blas.Trans or blas.ConjTrans

where a is an m×n band matrix kL subdiagonals and kU super-diagonals, and m and n refer to the size of the full dense matrix it represents. 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 computes

C = beta * C + alpha * A * B,

where A, B, and C are dense matrices, 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^T * 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 is a scalar.

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^T

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, 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 float64, b float64) (c float64, s float64, r float64, z float64)

func (Implementation) Drotm

func (Implementation) Drotm(n int, x []float64, incX int, y []float64, incY int, p blas.DrotmParams)

func (Implementation) Drotmg

func (Implementation) Drotmg(d1 float64, d2 float64, b1 float64, b2 float64) (p blas.DrotmParams, rd1 float64, rd2 float64, rb1 float64)

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

y = alpha * A * x + beta * y

where A is an n×n symmetric banded matrix, 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

func (Implementation) Dsdot(n int, x []float32, incX int, y []float32, incY int) float64

Dsdot computes the dot product of the two vectors

\sum_i x[i]*y[i]

func (Implementation) Dspmv

func (Implementation) Dspmv(ul blas.Uplo, n int, alpha float64, ap, x []float64, incX int, beta float64, y []float64, incY int)

Dspmv performs

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

func (Implementation) Dspr(ul blas.Uplo, n int, alpha float64, x []float64, incX int, ap []float64)

Dspr computes the rank-one operation

a += alpha * x * x^T

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^T + alpha * y * x^T,

where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha is a scalar.

func (Implementation) Dswap

func (Implementation) Dswap(n int, x []float64, incX int, y []float64, incY int)

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

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 computes

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 rank-one update

a += alpha * x * x^T

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^T + alpha * y * x^T

where A is a symmetric n×n matrix, x and y are vectors, and alpha is a scalar.

func (Implementation) Dsyr2k

func (Implementation) Dsyr2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha float64, a []float64, lda int, b []float64, ldb int, beta float64, c []float64, ldc int)

Dsyr2k performs the symmetric rank 2k operation

C = alpha * A * B^T + alpha * B * A^T + beta * C

where C is an n×n symmetric matrix. A and B are n×k matrices if tA == NoTrans and k×n otherwise. alpha and beta are scalars.

func (Implementation) Dsyrk

func (Implementation) Dsyrk(ul blas.Uplo, t blas.Transpose, n, k int, alpha float64, a []float64, lda int, beta float64, c []float64, ldc int)

Dsyrk performs the symmetric rank-k operation

C = alpha * A * A^T + beta*C

C is an n×n symmetric matrix. A is an n×k matrix if tA == blas.NoTrans, and a k×n matrix otherwise. 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 computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

where A is an n×n triangular banded matrix with k 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

A * x = b

where A is an n×n triangular banded matrix with k diagonals in packed format, and x is a vector. 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, x []float64, incX int)

Dtpmv computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

where A is an n×n unit 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, x []float64, incX int)

Dtpsv solves

A * x = b if tA == blas.NoTrans
A^T * x = b if tA == blas.Trans or blas.ConjTrans

where A is an n×n triangular matrix in packed format and x is a vector. 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

B = alpha * A * B,   if tA == blas.NoTrans and side == blas.Left,
B = alpha * A^T * 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^T, if tA == blas.Trans or blas.ConjTrans, and side == blas.Right,

where A is an n×n or m×m triangular matrix, and B is an m×n matrix.

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 computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

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

A * X = alpha * B,   if tA == blas.NoTrans side == blas.Left,
A^T * X = alpha * B, if tA == blas.Trans or blas.ConjTrans, and side == blas.Left,
X * A = alpha * B,   if tA == blas.NoTrans side == blas.Right,
X * A^T = 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 is an m×n matrix, 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

A * x = b if tA == blas.NoTrans
A^T * x = b if tA == blas.Trans or blas.ConjTrans

A is an n×n triangular matrix and x is a vector. 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

func (Implementation) Dznrm2

func (Implementation) Dznrm2(n int, x []complex128, incX int) float64

func (Implementation) Icamax

func (Implementation) Icamax(n int, x []complex64, incX int) int

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.

func (Implementation) Izamax

func (Implementation) Izamax(n int, x []complex128, incX int) int

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.

func (Implementation) Saxpy

func (Implementation) Saxpy(n int, alpha float32, x []float32, incX int, y []float32, incY int)

Saxpy adds alpha times x to y

y[i] += alpha * x[i] for all i

func (Implementation) Scasum

func (Implementation) Scasum(n int, x []complex64, incX int) float32

func (Implementation) Scnrm2

func (Implementation) Scnrm2(n int, x []complex64, incX int) float32

func (Implementation) Scopy

func (Implementation) Scopy(n int, x []float32, incX int, y []float32, incY int)

Scopy copies the elements of x into the elements of y.

y[i] = x[i] for all i

func (Implementation) Sdot

func (Implementation) Sdot(n int, x []float32, incX int, y []float32, incY int) float32

Sdot computes the dot product of the two vectors

\sum_i x[i]*y[i]

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]

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 computes

y = alpha * A * x + beta * y if tA == blas.NoTrans
y = alpha * A^T * x + beta * y if tA == blas.Trans or blas.ConjTrans

where a is an m×n band matrix kL subdiagonals and kU super-diagonals, and m and n refer to the size of the full dense matrix it represents. x and y are vectors, and alpha and beta are scalars.

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 computes

C = beta * C + alpha * A * B,

where A, B, and C are dense matrices, and alpha and beta are scalars. tA and tB specify whether A or B are transposed.

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^T * 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 is a scalar.

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^T

where A is an m×n dense matrix, x and y are vectors, and alpha is a scalar.

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.

func (Implementation) Srot

func (Implementation) Srot(n int, x []float32, incX int, y []float32, incY int, c, s float32)

Srot applies a plane transformation.

x[i] = c * x[i] + s * y[i]
y[i] = c * y[i] - s * x[i]

func (Implementation) Srotg

func (Implementation) Srotg(a float32, b float32) (c float32, s float32, r float32, z float32)

func (Implementation) Srotm

func (Implementation) Srotm(n int, x []float32, incX int, y []float32, incY int, p blas.SrotmParams)

func (Implementation) Srotmg

func (Implementation) Srotmg(d1 float32, d2 float32, b1 float32, b2 float32) (p blas.SrotmParams, rd1 float32, rd2 float32, rb1 float32)

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

y = alpha * A * x + beta * y

where A is an n×n symmetric banded matrix, x and y are vectors, and alpha and beta are scalars.

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.

func (Implementation) Sspmv

func (Implementation) Sspmv(ul blas.Uplo, n int, alpha float32, ap, x []float32, incX int, beta float32, y []float32, incY int)

Sspmv performs

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

func (Implementation) Sspr(ul blas.Uplo, n int, alpha float32, x []float32, incX int, ap []float32)

Sspr computes the rank-one operation

a += alpha * x * x^T

where a is an n×n symmetric matrix in packed format, x is a vector, and alpha is a scalar.

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^T + alpha * y * x^T,

where A is an n×n symmetric matrix in packed format, x and y are vectors, and alpha is a scalar.

func (Implementation) Sswap

func (Implementation) Sswap(n int, x []float32, incX int, y []float32, incY int)

Sswap exchanges the elements of two vectors.

x[i], y[i] = y[i], x[i] for all i

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

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

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

func (Implementation) Ssyr(ul blas.Uplo, n int, alpha float32, x []float32, incX int, a []float32, lda int)

Ssyr performs the rank-one update

a += alpha * x * x^T

where a is an n×n symmetric matrix, and x is a vector.

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^T + alpha * y * x^T

where A is a symmetric n×n matrix, x and y are vectors, and alpha is a scalar.

func (Implementation) Ssyr2k

func (Implementation) Ssyr2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha float32, a []float32, lda int, b []float32, ldb int, beta float32, c []float32, ldc int)

Ssyr2k performs the symmetric rank 2k operation

C = alpha * A * B^T + alpha * B * A^T + beta * C

where C is an n×n symmetric matrix. A and B are n×k matrices if tA == NoTrans and k×n otherwise. alpha and beta are scalars.

func (Implementation) Ssyrk

func (Implementation) Ssyrk(ul blas.Uplo, t blas.Transpose, n, k int, alpha float32, a []float32, lda int, beta float32, c []float32, ldc int)

Ssyrk performs the symmetric rank-k operation

C = alpha * A * A^T + beta*C

C is an n×n symmetric matrix. A is an n×k matrix if tA == blas.NoTrans, and a k×n matrix otherwise. alpha and beta are scalars.

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 computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

where A is an n×n triangular banded matrix with k diagonals, and x is a vector.

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

A * x = b

where A is an n×n triangular banded matrix with k diagonals in packed format, and x is a vector. 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) Stpmv

func (Implementation) Stpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []float32, incX int)

Stpmv computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

where A is an n×n unit triangular matrix in packed format, and x is a vector.

func (Implementation) Stpsv

func (Implementation) Stpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []float32, incX int)

Stpsv solves

A * x = b if tA == blas.NoTrans
A^T * x = b if tA == blas.Trans or blas.ConjTrans

where A is an n×n triangular matrix in packed format and x is a vector. 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) 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

B = alpha * A * B,   if tA == blas.NoTrans and side == blas.Left,
B = alpha * A^T * 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^T, if tA == blas.Trans or blas.ConjTrans, and side == blas.Right,

where A is an n×n or m×m triangular matrix, and B is an m×n matrix.

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 computes

x = A * x if tA == blas.NoTrans
x = A^T * x if tA == blas.Trans or blas.ConjTrans

A is an n×n Triangular matrix and x is a vector.

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

A * X = alpha * B,   if tA == blas.NoTrans side == blas.Left,
A^T * X = alpha * B, if tA == blas.Trans or blas.ConjTrans, and side == blas.Left,
X * A = alpha * B,   if tA == blas.NoTrans side == blas.Right,
X * A^T = 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 is an m×n matrix, 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) 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

A * x = b if tA == blas.NoTrans
A^T * x = b if tA == blas.Trans or blas.ConjTrans

A is an n×n triangular matrix and x is a vector. 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) Zaxpy

func (Implementation) Zaxpy(n int, alpha complex128, x []complex128, incX int, y []complex128, incY int)

func (Implementation) Zcopy

func (Implementation) Zcopy(n int, x []complex128, incX int, y []complex128, incY int)

func (Implementation) Zdotc

func (Implementation) Zdotc(n int, x []complex128, incX int, y []complex128, incY int) (dotc complex128)

func (Implementation) Zdotu

func (Implementation) Zdotu(n int, x []complex128, incX int, y []complex128, incY int) (dotu complex128)

func (Implementation) Zdscal

func (Implementation) Zdscal(n int, alpha float64, x []complex128, incX int)

func (Implementation) Zgbmv

func (Implementation) Zgbmv(tA blas.Transpose, m, n, kL, kU int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)

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)

func (Implementation) Zgemv

func (Implementation) Zgemv(tA blas.Transpose, m, n int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)

func (Implementation) Zgerc

func (Implementation) Zgerc(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)

func (Implementation) Zgeru

func (Implementation) Zgeru(m, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)

func (Implementation) Zhbmv

func (Implementation) Zhbmv(ul blas.Uplo, n, k int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)

func (Implementation) Zhemm

func (Implementation) Zhemm(s blas.Side, ul blas.Uplo, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)

func (Implementation) Zhemv

func (Implementation) Zhemv(ul blas.Uplo, n int, alpha complex128, a []complex128, lda int, x []complex128, incX int, beta complex128, y []complex128, incY int)

func (Implementation) Zher

func (Implementation) Zher(ul blas.Uplo, n int, alpha float64, x []complex128, incX int, a []complex128, lda int)

func (Implementation) Zher2

func (Implementation) Zher2(ul blas.Uplo, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, a []complex128, lda int)

func (Implementation) Zher2k

func (Implementation) Zher2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta float64, c []complex128, ldc int)

func (Implementation) Zherk

func (Implementation) Zherk(ul blas.Uplo, t blas.Transpose, n, k int, alpha float64, a []complex128, lda int, beta float64, c []complex128, ldc int)

func (Implementation) Zhpmv

func (Implementation) Zhpmv(ul blas.Uplo, n int, alpha complex128, ap, x []complex128, incX int, beta complex128, y []complex128, incY int)

func (Implementation) Zhpr

func (Implementation) Zhpr(ul blas.Uplo, n int, alpha float64, x []complex128, incX int, ap []complex128)

func (Implementation) Zhpr2

func (Implementation) Zhpr2(ul blas.Uplo, n int, alpha complex128, x []complex128, incX int, y []complex128, incY int, ap []complex128)

func (Implementation) Zscal

func (Implementation) Zscal(n int, alpha complex128, x []complex128, incX int)

func (Implementation) Zswap

func (Implementation) Zswap(n int, x []complex128, incX int, y []complex128, incY int)

func (Implementation) Zsymm

func (Implementation) Zsymm(s blas.Side, ul blas.Uplo, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)

func (Implementation) Zsyr2k

func (Implementation) Zsyr2k(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, b []complex128, ldb int, beta complex128, c []complex128, ldc int)

func (Implementation) Zsyrk

func (Implementation) Zsyrk(ul blas.Uplo, t blas.Transpose, n, k int, alpha complex128, a []complex128, lda int, beta complex128, c []complex128, ldc int)

func (Implementation) Ztbmv

func (Implementation) Ztbmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []complex128, lda int, x []complex128, incX int)

func (Implementation) Ztbsv

func (Implementation) Ztbsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n, k int, a []complex128, lda int, x []complex128, incX int)

func (Implementation) Ztpmv

func (Implementation) Ztpmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []complex128, incX int)

func (Implementation) Ztpsv

func (Implementation) Ztpsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, ap, x []complex128, incX int)

func (Implementation) Ztrmm

func (Implementation) Ztrmm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int)

func (Implementation) Ztrmv

func (Implementation) Ztrmv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []complex128, lda int, x []complex128, incX int)

func (Implementation) Ztrsm

func (Implementation) Ztrsm(s blas.Side, ul blas.Uplo, tA blas.Transpose, d blas.Diag, m, n int, alpha complex128, a []complex128, lda int, b []complex128, ldb int)

func (Implementation) Ztrsv

func (Implementation) Ztrsv(ul blas.Uplo, tA blas.Transpose, d blas.Diag, n int, a []complex128, lda int, x []complex128, incX int)

Notes

Bugs

  • The cgo package is intrinsically dependent on the underlying C implementation. The BLAS standard is silent on a number of behaviors, including but not limited to how NaN values are treated. For this reason the result of computations performed by the cgo BLAS package may disagree with the results produced by the native BLAS package. The cgo package is tested against OpenBLAS; use of other backing BLAS C libraries may result in test failure because of this.

Jump to

Keyboard shortcuts

? : This menu
/ : Search site
f or F : Jump to
y or Y : Canonical URL