blas

package
v0.0.0-...-3dbd70c Latest Latest
Warning

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

Go to latest
Published: Nov 17, 2014 License: GPL-3.0 Imports: 5 Imported by: 0

Documentation

Overview

package blas provides an interface to BLAS (Basic Linear Algebra Subprograms) provided by the GSL (GNU Scientific Library).

Matrix representation

Matrices are represented by 2-dimensional slices like [][]float32. However, BLAS requires contiguous matrix memory. Therefore, matrices must not be manually allocated, one has to use:

MakeFloat32Matrix(rows, cols)
MakeFloat64Matrix(rows, cols)
MakeComplex64Matrix(rows, cols)
MakeComplex128Matrix(rows, cols)

The elements of these matrices can be manipulated as usual. However, the slice values themselves (length, capacity, reference to data) should not be changed. E.g.: rows must not be swapped with constructs like M[0], M[1] = M[1], M[0]. Blas functions may panic or return unexpected results when passed improperly constructed matrices, but memory safety is in principle guaranteed.

Submatrices can be constructed with the *Submatrix functions. These are also safe to pass to BLAS functions. They share storage with the original matrix.

Just like the GLS BLAS interface this package was derived from, operations on packed/band matrices are not supported. These are available through the low-level cblas package.

The library assumes that arrays, vectors and matrices passed as modifiable arguments are not aliased and do not overlap with each other. This removes the need for the library to handle overlapping memory regions as a special case, and allows additional optimizations to be used. If overlapping memory regions are passed as modifiable arguments then the results of such functions will be undefined. If the arguments will not be modified then overlapping or aliased memory regions can be safely used.

Naming

Each routine has a name which specifies the operation, the type of matrices involved and their precisions. Some of the most common operations and their names are given below,

DOT  scalar product, x^T y
AXPY vector sum, \alpha x + y
MV   matrix-vector product, A x
SV   matrix-vector solve, inv(A) x
MM   matrix-matrix product, A B
SM   matrix-matrix solve, inv(A) B

The types of matrices are,

GE general
GB general band
SY symmetric
SB symmetric band
SP symmetric packed
HE hermitian
HB hermitian band
HP hermitian packed
TR triangular
TB triangular band
TP triangular packed

Each operation is defined for four precisions,

S single real
D double real
C single complex
Z double complex

Thus, for example, the name SGEMM stands for “single-precision general matrix-matrix multiply” and ZGEMM stands for “double-precision complex matrix-matrix multiply”.

Index

Examples

Constants

This section is empty.

Variables

This section is empty.

Functions

func CAXPY

func CAXPY(alpha complex64, X []complex64, incX int, Y []complex64, incY int)

Replaces Y by (alpha*X) + Y. Every incX'th and incY'th element is used.

Example
alpha := complex64(complex(0, 1))
X := []complex64{1, 2, 3}
incX := 1
Y := []complex64{4, 5, 6}
incY := 1
CAXPY(alpha, X, incX, Y, incY)
fmt.Println(Y)
Output:

[(4+1i) (5+2i) (6+3i)]

func CCOPY

func CCOPY(X []complex64, incX int, Y []complex64, incY int)

Copies X to Y. Every incX'th and incY'th element is used.

Example
X := []complex64{1, 666, 2, 666}
incX := 2
Y := []complex64{-1, -2}
incY := 1
CCOPY(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[(1+0i) (666+0i) (2+0i) (666+0i)] [(1+0i) (2+0i)]

func CDOTC

func CDOTC(X []complex64, incX int, Y []complex64, incY int) complex64

Calculates the dot product of the complex conjugate of X with Y. Every incX'th and incY'th element is used.

Example
X := []complex64{complex(1, 0), complex(1, 1)}
incX := 1
Y := []complex64{complex(0, 0), complex(2, 0)}
incY := 1
result := CDOTC(X, incX, Y, incY)
fmt.Println(result)
Output:

(2-2i)

func CDOTU

func CDOTU(X []complex64, incX int, Y []complex64, incY int) complex64

Computes the dot product of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []complex64{complex(1, 0), complex(1, 1)}
incX := 1
Y := []complex64{complex(0, 0), complex(2, 0)}
incY := 1
result := CDOTU(X, incX, Y, incY)
fmt.Println(result)
Output:

(2+2i)

func CGEMM

func CGEMM(transA Transpose, transB Transpose, alpha complex64, A [][]complex64, B [][]complex64, beta complex64, C [][]complex64)

General matrix-matrix multiplication:

C = alpha*A*B + beta*C

where A and B can optionally be transposed by specifying transA, transB = Trans/NoTrans

Example
alpha := complex64(1.0)
A := matrix.MakeComplex64(2, 4)
A[0][0] = 1
A[0][1] = complex(0, -1)
A[1][1] = 2

B := matrix.MakeComplex64(4, 3)
B[1][0] = 2
B[2][1] = 3

beta := complex64(0.0)
C := matrix.MakeComplex64(2, 3)
CGEMM(NoTrans, NoTrans, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[(1+0i) (0-1i) (0+0i) (0+0i)] [(0+0i) (2+0i) (0+0i) (0+0i)]] * [[(0+0i) (0+0i) (0+0i)] [(2+0i) (0+0i) (0+0i)] [(0+0i) (3+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]] = [[(0-2i) (0+0i) (0+0i)] [(4+0i) (0+0i) (0+0i)]]

func CGEMV

func CGEMV(transA Transpose, alpha complex64, A [][]complex64, X []complex64, incX int, beta complex64, Y []complex64, incY int)

Matrix-vector multiplication plus vector with optional matrix transpose. When transA == NoTrans this computes:

Y = alpha*A*X + beta*Y

When transA == Trans this computes:

Y = alpha*(A^T)*X + beta*Y

Every incX'th element of X and incY'th element of Y is used. Matrices must be allocated with MakeComplex64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeComplex64(2, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []complex64{-1, 4, 0}
incX := 1
Y := []complex64{0, 0}
incY := 1

alpha := complex64(complex(1, 0))
beta := complex64(complex(0, 0))

CGEMV(NoTrans, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)

CGEMV(Trans, alpha, A, Y, incX, beta, X, incY)
fmt.Println(A, "^T*", Y, "=", X)
Output:

[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)]] * [(-1+0i) (4+0i) (0+0i)] = [(4+0i) (-3+0i)]
[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)]] ^T* [(4+0i) (-3+0i)] = [(-9+0i) (4+0i) (8+0i)]

func CGERC

func CGERC(alpha complex64, X []complex64, incX int, Y []complex64, incY int, A [][]complex64)

Computes the conjugate rank-1 update:

A = alpha*X*conj(Y^T) + A

Matrices must be allocated with MakeComplex64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := complex64(1)
X := []complex64{1, 2}
incX := 1
Y := []complex64{complex(0, 3), complex(0, 4), complex(0, 5)}
incY := 1
A := matrix.MakeComplex64(2, 3)
CGERC(alpha, X, incX, Y, incY, A)
fmt.Println(X, "* conj", Y, "^T = ", A)
Output:

[(1+0i) (2+0i)] * conj [(0+3i) (0+4i) (0+5i)] ^T =  [[(0-3i) (0-4i) (0-5i)] [(0-6i) (0-8i) (0-10i)]]

func CGERU

func CGERU(alpha complex64, X []complex64, incX int, Y []complex64, incY int, A [][]complex64)

Computes the rank-1 update:

A = alpha*X*(Y^T) + A

Matrices must be allocated with MakeComplex64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := complex64(1)
X := []complex64{1, 2}
incX := 1
Y := []complex64{3, 4, 5}
incY := 1
A := matrix.MakeComplex64(2, 3)
CGERU(alpha, X, incX, Y, incY, A)
fmt.Println(X, "*", Y, "^T = ", A)
Output:

[(1+0i) (2+0i)] * [(3+0i) (4+0i) (5+0i)] ^T =  [[(3+0i) (4+0i) (5+0i)] [(6+0i) (8+0i) (10+0i)]]

func CHEMM

func CHEMM(side Side, uplo Uplo, alpha complex64, A [][]complex64, B [][]complex64, beta complex64, C [][]complex64)

Computes the matrix-matrix product and sum

C = alpha A B + beta C (for Side==Left)
C = alpha B A + beta C (for Side==Right)

where the matrix A is hermitian. When Uplo is Upper then the upper triangle and diagonal of A are used, and when Uplo is Lower then the lower triangle and diagonal of A are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := complex64(1.0)
A := matrix.MakeComplex64(2, 2)
A[0][0] = 1
A[0][1] = complex(0, -1)
A[1][1] = 2

B := matrix.MakeComplex64(2, 2)
B[1][0] = 2
B[0][1] = 3

beta := complex64(0.0)
C := matrix.MakeComplex64(2, 2)
CHEMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(0-2i) (3+0i)] [(4+0i) (0+3i)]]

func CHEMV

func CHEMV(uplo Uplo, alpha complex64, A [][]complex64, X []complex64, incX int, beta complex64, Y []complex64, incY int)

Hermitian matrix-vector product:

Y = alpha*X + beta*Y

Matrices must be allocated with MakeComplex64Matrix to ensure contiguous underlying storage, With uplo == Upper or Lower, the upper or lower triangular part of A is used. Every incX'th and incY'th element is used.

Example
A := matrix.MakeComplex64(3, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []complex64{-1, 4, 0}
incX := 1
Y := []complex64{0, 0, 0}
incY := 1

alpha := complex64(complex(1, 0))
beta := complex64(complex(0, 0))

CHEMV(Upper, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A)
Output:

[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func CHER

func CHER(uplo Uplo, alpha float32, X []complex64, incX int, A [][]complex64)

Computes the hermitian rank-1 update:

A = alpha*X*conj(X^T) + A

A is a hermitian matrix. Only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := float32(1.0)
X := []complex64{complex(0, 1), 2}
incX := 1
A := matrix.MakeComplex64(2, 2)
CHER(Upper, alpha, X, incX, A)
fmt.Println(A)
Output:

[[(1+0i) (0+2i)] [(0+0i) (4+0i)]]

func CHER2

func CHER2(uplo Uplo, alpha complex64, X []complex64, incX int, Y []complex64, incY int, A [][]complex64)

Computes the hermitian rank-2 update: A = alpha X conj(Y^T) + conj(alpha) Y conj(X^T) + A A is a hermitian matrix. Only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex64(1.0)
X := []complex64{complex(0, 1), 2}
incX := 1
Y := []complex64{2, complex(3, 1)}
incY := 1
A := matrix.MakeComplex64(2, 2)
CHER2(Upper, alpha, X, incX, Y, incY, A)
fmt.Println(A)
Output:

[[(2+0i) (0+4i)] [(0+0i) (8+0i)]]

func CHER2K

func CHER2K(uplo Uplo, trans Transpose, alpha complex64, A [][]complex64, B [][]complex64, beta float32, C [][]complex64)

Computes a rank-2k update of the hermitian matrix C,

C = alpha A B^H + alpha^* B A^H + beta C (for trans==NoTrans)
C = alpha A^H B + alpha^* B^H A + beta C (for Trans==ConjTrans)

Since the matrix C is hermitian only its upper half or lower half need to be stored. When Uplo is Upper then the upper triangle and diagonal of C are used, and when Uplo is Lower then the lower triangle and diagonal of C are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := complex64(complex(2.0, 3.0))
A := matrix.MakeComplex64(2, 2)
A[1][1] = complex(0, 1)
B := matrix.MakeComplex64(2, 2)
B[0][1] = complex(1, 0)
beta := float32(2.0)
C := matrix.MakeComplex64(2, 2)
CHER2K(Upper, NoTrans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(0+0i) (-3-2i)] [(0+0i) (0+0i)]]

func CHERK

func CHERK(uplo Uplo, trans Transpose, alpha float32, A [][]complex64, beta float32, C [][]complex64)

Computes a rank-k update of the hermitian matrix C

C = alpha A A^H + beta C (for trans==NoTrans)
C = alpha A^H A + beta C (for trans==ConjTrans)

Since the matrix C is hermitian only its upper half or lower half need to be stored. When Uplo is Upper then the upper triangle and diagonal of C are used, and when Uplo is Lower then the lower triangle and diagonal of C are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := float32(1.0)
A := matrix.MakeComplex64(2, 2)
A[1][1] = complex(0, 1)
beta := float32(0.0)
C := matrix.MakeComplex64(2, 2)
CHERK(Upper, ConjTrans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[(0+0i) (0+0i)] [(0+0i) (1+0i)]]

func CSCAL

func CSCAL(alpha complex64, X []complex64, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := complex64(complex(0, 1))
X := []complex64{1, 666, 2, 666}
incX := 2
CSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[(0+1i) (666+0i) (0+2i) (666+0i)]

func CSSCAL

func CSSCAL(alpha float32, X []complex64, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := float32(2)
X := []complex64{1, 666, 2, 666}
incX := 2
CSSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[(2+0i) (666+0i) (4+0i) (666+0i)]

func CSWAP

func CSWAP(X []complex64, incX int, Y []complex64, incY int)

Exchanges the elements of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []complex64{1, 666, 2, 666}
incX := 2
Y := []complex64{-1, -2}
incY := 1
CSWAP(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[(-1+0i) (666+0i) (-2+0i) (666+0i)] [(1+0i) (2+0i)]

func CSYMM

func CSYMM(side Side, uplo Uplo, alpha complex64, A [][]complex64, B [][]complex64, beta complex64, C [][]complex64)

Computes the matrix-matrix product

C = alpha*A*B + beta*C  (for Side==Left)
C = alpha*B*A + beta*C  (for Side==Right)

where the matrix A is symmetric. Only its upper half or lower half is used, specified by uplo=Upper/Lower.

Example
alpha := complex64(1)
A := matrix.MakeComplex64(3, 3)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeComplex64(3, 2)
B[1][0] = 2
B[2][1] = 3

beta := complex64(0)
C := matrix.MakeComplex64(3, 2)
CSYMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[(1+0i) (-1+0i) (0+0i)] [(0+0i) (2+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]] * [[(0+0i) (0+0i)] [(2+0i) (0+0i)] [(0+0i) (3+0i)]] = [[(-2+0i) (0+0i)] [(4+0i) (0+0i)] [(0+0i) (0+0i)]]

func CSYR2K

func CSYR2K(uplo Uplo, trans Transpose, alpha complex64, A [][]complex64, B [][]complex64, beta complex64, C [][]complex64)

Computes a rank-2k update of the symmetric matrix C:

C = alpha A*(B^T) + alpha B*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*B + alpha (B^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex64(1.0)
A := matrix.MakeComplex64(2, 3)
A[0][0] = 1

B := matrix.MakeComplex64(2, 3)
B[0][0] = 1

beta := complex64(0.0)
C := matrix.MakeComplex64(3, 3)

CSYR2K(Upper, Trans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(2+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func CSYRK

func CSYRK(uplo Uplo, trans Transpose, alpha complex64, A [][]complex64, beta complex64, C [][]complex64)

Computes a rank-k update of the symmetric matrix C:

C = alpha A*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex64(1)
A := matrix.MakeComplex64(2, 3)
A[0][0] = 1

beta := complex64(0)
C := matrix.MakeComplex64(3, 3)

CSYRK(Upper, Trans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[(1+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func CTRMM

func CTRMM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha complex64, A [][]complex64, B [][]complex64)

Computes the matrix-matrix product

B = alpha op(A) B (for Side is Left)
B = alpha B op(A) (for Side is Right)

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := complex64(1)
A := matrix.MakeComplex64(3, 3)
A[1][1] = 1
A[1][2] = complex(0, 1)
B := matrix.MakeComplex64(2, 3)
B[1][2] = 1

CTRMM(Right, Upper, ConjTrans, NonUnit, alpha, A, B)
fmt.Println(B)
Output:

[[(0+0i) (0+0i) (0+0i)] [(0+0i) (0-1i) (0+0i)]]

func CTRMV

func CTRMV(uplo Uplo, transA Transpose, diag Diag, A [][]complex64, X []complex64, incX int)

Triangular matrix-vector multiplication plus vector with optional matrix transpose. With uplo == Upper or Lower, the upper or lower triangular part of A is used. diag specifies whether the matrix is unit triangular (). When transA == NoTrans this computes:

X = A*X

When transA == Trans this computes:

X = (A^T)*X

Every incX'th element of X is used. Matrices must be allocated with MakeComplex64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeComplex64(2, 2)
A[0][0] = 1
A[1][1] = 1

A[0][1] = 2

X := []complex64{-1, 4}
incX := 1

fmt.Println(A, "*", X)
CTRMV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println("=", X)
Output:

[[(1+0i) (2+0i)] [(0+0i) (1+0i)]] * [(-1+0i) (4+0i)]
= [(7+0i) (4+0i)]

func CTRSM

func CTRSM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha complex64, A [][]complex64, B [][]complex64)

Computes the inverse-matrix matrix product

B = alpha op(inv(A))B for Side is Left
B = alpha B op(inv(A)) for Side is Right.

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := complex64(1)
A := matrix.MakeComplex64(3, 3)
A[1][1] = 1
A[1][2] = complex(0, 1)
B := matrix.MakeComplex64(2, 3)
B[1][2] = 1

CTRSM(Right, Upper, ConjTrans, Unit, alpha, A, B)
fmt.Println(B)
Output:

[[(0+0i) (0+0i) (0+0i)] [(0+0i) (0+1i) (1+0i)]]

func CTRSV

func CTRSV(uplo Uplo, transA Transpose, diag Diag, A [][]complex64, X []complex64, incX int)

Computes

inv(op(A)) x for x

where op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of the matrix is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
// solve:
//  x - y = 2
//      y = 3
// for x, y:
A := matrix.MakeComplex64(2, 2)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 1
X := []complex64{2, 3}
incX := 1
CTRSV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println(X)
Output:

[(5+0i) (3+0i)]

func DASUM

func DASUM(X []float64, incX int) float64

Computes the sum of the absolute values of the elements in vector X. Every incX'th element is used.

Example
X := []float64{1, 666, 2, 666}
incX := 2
result := DASUM(X, incX)
fmt.Println(result)
Output:

3

func DAXPY

func DAXPY(alpha float64, X []float64, incX int, Y []float64, incY int)

Replaces Y by (alpha*X) + Y. Every incX'th and incY'th element is used.

Example
alpha := float64(2)
X := []float64{1, 2, 3}
incX := 1
Y := []float64{4, 5, 6}
incY := 1
DAXPY(alpha, X, incX, Y, incY)
fmt.Println(Y)
Output:

[6 9 12]

func DCOPY

func DCOPY(X []float64, incX int, Y []float64, incY int)

Copies X to Y. Every incX'th and incY'th element is used.

Example
X := []float64{1, 666, 2, 666}
incX := 2
Y := []float64{-1, -2}
incY := 1
DCOPY(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[1 666 2 666] [1 2]

func DDOT

func DDOT(X []float64, incX int, Y []float64, incY int) float64

Computes the dot product of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []float64{2, 0}
incX := 1
Y := []float64{3, 5}
incY := 1
result := DDOT(X, incX, Y, incY)
fmt.Println(result)
Output:

6

func DGEMM

func DGEMM(transA Transpose, transB Transpose, alpha float64, A [][]float64, B [][]float64, beta float64, C [][]float64)

General matrix-matrix multiplication:

C = alpha*A*B + beta*C

where A and B can optionally be transposed by specifying transA, transB = Trans/NoTrans

Example
alpha := 1.0
A := matrix.MakeFloat64(2, 4)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeFloat64(4, 3)
B[1][0] = 2
B[2][1] = 3

beta := 0.0
C := matrix.MakeFloat64(2, 3)
DGEMM(NoTrans, NoTrans, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[1 -1 0 0] [0 2 0 0]] * [[0 0 0] [2 0 0] [0 3 0] [0 0 0]] = [[-2 0 0] [4 0 0]]

func DGEMV

func DGEMV(transA Transpose, alpha float64, A [][]float64, X []float64, incX int, beta float64, Y []float64, incY int)

Matrix-vector multiplication plus vector with optional matrix transpose. When transA == NoTrans this computes:

Y = alpha*A*X + beta*Y

When transA == Trans this computes:

Y = alpha*(A^T)*X + beta*Y

Every incX'th element of X and incY'th element of Y is used. Matrices must be allocated with MakeFloat64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeFloat64(2, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []float64{-1, 4, 0}
incX := 1
Y := []float64{0, 0}
incY := 1

alpha := 1.0
beta := 0.0

DGEMV(NoTrans, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)

DGEMV(Trans, alpha, A, Y, incX, beta, X, incY)
fmt.Println(A, "^T*", Y, "=", X)
Output:

[[0 1 2] [3 0 0]] * [-1 4 0] = [4 -3]
[[0 1 2] [3 0 0]] ^T* [4 -3] = [-9 4 8]

func DGER

func DGER(alpha float64, X []float64, incX int, Y []float64, incY int, A [][]float64)

Computes

A = alpha*X*(Y^T) + A

Matrices must be allocated with MakeFloat64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := 1.0
X := []float64{1, 2}
incX := 1
Y := []float64{3, 4, 5}
incY := 1
A := matrix.MakeFloat64(2, 3)
DGER(alpha, X, incX, Y, incY, A)
fmt.Println(X, "*", Y, "^T = ", A)
Output:

[1 2] * [3 4 5] ^T =  [[3 4 5] [6 8 10]]

func DNRM2

func DNRM2(X []float64, incX int) float64

Computes the L2 norm (Euclidian length) of vector X. Every incX'th element is used.

Example
X := []float64{1, 666, 1, 666}
incX := 2
result := DNRM2(X, incX)
fmt.Println(result)
Output:

1.4142135623730951

func DROT

func DROT(X []float64, incX int, Y []float64, incY int, c float64, s float64)

Applies a plane rotation to the two vectors X and Y. This routine computes:

[ x_i ]   [ c s ] [ x_i ]
[ y_i ] = [-s c ] [ y_i ]

for all i with strides incX, incY.

Example
X := []float64{1, 0, 1}
incX := 1
Y := []float64{0, 1, 1}
incY := 1
theta := math.Pi / 4
c := math.Cos(theta)
s := math.Sin(theta)
DROT(X, incX, Y, incY, c, s)
fmt.Println(X, Y)
Output:

[0.7071067811865476 0.7071067811865475 1.414213562373095] [-0.7071067811865475 0.7071067811865476 1.1102230246251565e-16]

func DROTG

func DROTG(a float64, b float64) (r, c, s float64)

Constructs a Givens rotation matrix. Computes the values of c and s so that:

[ c  s ] [ a ]    [ r ]
[ -s c ] [ b ] =  [ 0 ].

Returns r, c and s.

func DROTM

func DROTM(X []float64, incX int, Y []float64, incY int, P [5]float64)

Applies the modified-Givens rotation of the vectors X and Y:

[ x_i ]   [ h_11 h_12 ] [ x_i ]
[ y_i ] = [ h_21 h_22 ] [ y_i ]

for all i with strides incX, incY

H stored in P according to the encoding used by DROTMG.

Example
d1 := 1.0
d2 := 1.0
b1 := 1.0
b2 := 1.0
P := DROTMG(d1, d2, b1, b2)

X := []float64{1}
incX := 1
Y := []float64{1}
incY := 1
DROTM(X, incX, Y, incY, P)
fmt.Println(X, Y)
Output:

[2] [0]

func DROTMG

func DROTMG(d1 float64, d2 float64, b1 float64, b2 float64) [5]float64

Constructs a modified Givens rotation matrix. The input scalars d1, d2, x1 and y1 define a 2-vector [a1 a2]' such that

[ b1 ]   [ d1^{1/2}  0      ] [ x1 ]
[ b2 ] = [   0     d2^{1/2} ] [ y1 ].

This subroutine determines the modified Givens rotation matrix H that transforms y1 and thus a2 to zero. A representation of this matrix is returned:

P[0]: defines the form of matrix H:

-2.0: matrix H contains the identity matrix.
-1.0: matrix H is identical to matrix SH (defined by the remaining values in the vector).
 0.0: H[1,2] and H[2,1] are obtained from matrix SH. The remaining values are both 1.0.
 1.0: H[1,1] and H[2,2] are obtained from matrix SH. H[1,2] is 1.0. H[2,1] is -1.0.

The other elements contain SH:

P[1]	contains SH[1,1].
P[2]	contains SH[2,1].
P[3]	contains SH[1,2].
P[4]	contains SH[2,2].
Example
d1 := 1.0
d2 := 1.0
b1 := 1.0
b2 := 1.0
result := DROTMG(d1, d2, b1, b2)
fmt.Println(result)
Output:

[1 1 0 0 1]

func DSCAL

func DSCAL(alpha float64, X []float64, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := 2.0
X := []float64{1, 666, 2, 666}
incX := 2
DSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[2 666 4 666]

func DSDOT

func DSDOT(X []float32, incX int, Y []float32, incY int) float64

Computes the dot product of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []float32{2, 0}
incX := 1
Y := []float32{3, 5}
incY := 1
result := DSDOT(X, incX, Y, incY)
fmt.Println(result)
Output:

6

func DSWAP

func DSWAP(X []float64, incX int, Y []float64, incY int)

Exchanges the elements of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []float64{1, 666, 2, 666}
incX := 2
Y := []float64{-1, -2}
incY := 1
DSWAP(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[-1 666 -2 666] [1 2]

func DSYMM

func DSYMM(side Side, uplo Uplo, alpha float64, A [][]float64, B [][]float64, beta float64, C [][]float64)

Computes the matrix-matrix product

C = alpha*A*B + beta*C  (for Side==Left)
C = alpha*B*A + beta*C  (for Side==Right)

where the matrix A is symmetric. Only its upper half or lower half is used, specified by uplo=Upper/Lower.

Example
alpha := 1.0
A := matrix.MakeFloat64(3, 3)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeFloat64(3, 2)
B[1][0] = 2
B[2][1] = 3

beta := 0.0
C := matrix.MakeFloat64(3, 2)
DSYMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[1 -1 0] [0 2 0] [0 0 0]] * [[0 0] [2 0] [0 3]] = [[-2 0] [4 0] [0 0]]

func DSYMV

func DSYMV(uplo Uplo, alpha float64, A [][]float64, X []float64, incX int, beta float64, Y []float64, incY int)

Symmetric matrix-vector multiplication:

Y = alpha*A*X + beta*Y

A is a symmetric matrix, where uplo (Upper or Lower) determines which part of A is used. Every incX'th element of X and incY'th element of Y is used.

Example
A := matrix.MakeFloat64(2, 2)
A[0][0] = 1
A[1][1] = 2
A[0][1] = 3
A[1][0] = 3

X := []float64{-1, 4}
incX := 1
Y := []float64{0, 0}
incY := 1

alpha := float64(1.0)
beta := float64(0.0)

DSYMV(Upper, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)
Output:

[[1 3] [3 2]] * [-1 4] = [11 5]

func DSYR

func DSYR(uplo Uplo, alpha float64, X []float64, incX int, A [][]float64)

Computes

A = A + alpha*X*(X^T)

A is a symmetric matrix, where uplo (Upper or Lower) determines which part of A is used.

Example
alpha := 1.0
X := []float64{1, 2}
incX := 1
A := matrix.MakeFloat64(2, 2)
DSYR(Upper, alpha, X, incX, A)
fmt.Println(X, "*", X, "^T = ", A)
Output:

[1 2] * [1 2] ^T =  [[1 2] [0 4]]

func DSYR2

func DSYR2(uplo Uplo, alpha float64, X []float64, incX int, Y []float64, incY int, A [][]float64)

Computes the symmetric rank-2 update

A = \alpha x y^T + \alpha y x^T + A

of the symmetric matrix A. Since the matrix A is symmetric only its upper half or lower half need to be stored. When Uplo is Upper then the upper triangle and diagonal of A are used, and when Uplo is Lower then the lower triangle and diagonal of A are used.

Example
alpha := float64(1)
X := []float64{1, 2}
incX := 1
Y := []float64{2, 4}
incY := 1
A := matrix.MakeFloat64(2, 2)
DSYR2(Upper, alpha, X, incX, Y, incY, A)
fmt.Println(A)
Output:

[[4 8] [0 16]]

func DSYR2K

func DSYR2K(uplo Uplo, trans Transpose, alpha float64, A [][]float64, B [][]float64, beta float64, C [][]float64)

Computes a rank-2k update of the symmetric matrix C:

C = alpha A*(B^T) + alpha B*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*B + alpha (B^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := float64(1.0)
A := matrix.MakeFloat64(2, 3)
A[0][0] = 1

B := matrix.MakeFloat64(2, 3)
B[0][0] = 1

beta := float64(0.0)
C := matrix.MakeFloat64(3, 3)

DSYR2K(Upper, Trans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[2 0 0] [0 0 0] [0 0 0]]

func DSYRK

func DSYRK(uplo Uplo, trans Transpose, alpha float64, A [][]float64, beta float64, C [][]float64)

Computes a rank-k update of the symmetric matrix C:

C = alpha A*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := 1.0
A := matrix.MakeFloat64(2, 3)
A[0][0] = 1

beta := 0.0
C := matrix.MakeFloat64(3, 3)

DSYRK(Upper, Trans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[1 0 0] [0 0 0] [0 0 0]]

func DTRMM

func DTRMM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha float64, A [][]float64, B [][]float64)

Computes the matrix-matrix product

B = alpha op(A) B (for Side is Left)
B = alpha B op(A) (for Side is Right)

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := float64(1.0)
A := matrix.MakeFloat64(3, 3)
A[1][1] = 1
A[1][2] = 1
B := matrix.MakeFloat64(2, 3)
B[1][1] = 1

DTRMM(Right, Lower, Trans, NonUnit, alpha, A, B)
fmt.Println(B)
Output:

[[0 0 0] [0 1 0]]

func DTRMV

func DTRMV(uplo Uplo, transA Transpose, diag Diag, A [][]float64, X []float64, incX int)

Triangular matrix-vector multiplication plus vector with optional matrix transpose. With uplo == Upper or Lower, the upper or lower triangular part of A is used. diag specifies whether the matrix is unit triangular (). When transA == NoTrans this computes:

X = A*X

When transA == Trans this computes:

X = (A^T)*X

Every incX'th element of X is used. Matrices must be allocated with MakeFloat64Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeFloat64(2, 2)
A[0][0] = 1
A[1][1] = 1

A[0][1] = 2

X := []float64{-1, 4}
incX := 1

fmt.Println(A, "*", X)
DTRMV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println("=", X)
Output:

[[1 2] [0 1]] * [-1 4]
= [7 4]

func DTRSM

func DTRSM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha float64, A [][]float64, B [][]float64)

Computes the inverse-matrix matrix product

B = alpha op(inv(A))B for Side is Left
B = alpha B op(inv(A)) for Side is Right.

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := float64(1.0)
A := matrix.MakeFloat64(3, 3)
B := matrix.MakeFloat64(3, 2)
B[1][1] = 1

DTRSM(Left, Upper, NoTrans, Unit, alpha, A, B)
fmt.Println(B)
Output:

[[0 0] [0 1] [0 0]]

func DTRSV

func DTRSV(uplo Uplo, transA Transpose, diag Diag, A [][]float64, X []float64, incX int)

Computes

inv(op(A)) x for x

where op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of the matrix is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
// solve:
//  x - y = 2
//      y = 3
// for x, y:
A := matrix.MakeFloat64(2, 2)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 1
X := []float64{2, 3}
incX := 1
DTRSV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println(X)
Output:

[5 3]

func DZASUM

func DZASUM(X []complex128, incX int) float64

Computes the sum of the absolute values of real and imaginary parts of elements in vector X. Every incX'th element is used.

Example
X := []complex128{complex(1, 0), complex(1, 1)}
incX := 1
result := DZASUM(X, incX)
fmt.Println(result)
Output:

3

func DZNRM2

func DZNRM2(X []complex128, incX int) float64

Computes the unitary norm of vector X. Every incX'th element is used.

Example
X := []complex128{complex(1, 0), complex(0, 1)}
incX := 1
result := DZNRM2(X, incX)
fmt.Println(result)
Output:

1.4142135623730951

func ICAMAX

func ICAMAX(X []complex64, incX int) int

Returns the index of the element with the largest absolute value in vector X. Every incX'th element is used.

Example
X := []complex64{complex(2, 0), complex(2, 2), complex(0, 2)}
incX := 1
result := ICAMAX(X, incX)
fmt.Println(result)
Output:

1

func IDAMAX

func IDAMAX(X []float64, incX int) int

Returns the index of the element with the largest absolute value in vector X. Every incX'th element is used.

Example
X := []float32{1, -999, 2, 3, 999}
incX := 1
result := ISAMAX(X, incX)
fmt.Println(result)
Output:

1

func ISAMAX

func ISAMAX(X []float32, incX int) int

Returns the index of the element with the largest absolute value in vector X. Every incX'th element is used.

Example
X := []float32{1, -999, 2, 3, 999}
incX := 1
result := ISAMAX(X, incX)
fmt.Println(result)
Output:

1

func IZAMAX

func IZAMAX(X []complex128, incX int) int

Returns the index of the element with the largest absolute value in vector X. Every incX'th element is used.

Example
X := []complex128{complex(2, 0), complex(2, 2), complex(0, 2)}
incX := 1
result := IZAMAX(X, incX)
fmt.Println(result)
Output:

1

func SASUM

func SASUM(X []float32, incX int) float32

Computes the sum of the absolute values of the elements in vector X. Every incX'th element is used.

Example
X := []float32{1, 666, 2, 666}
incX := 2
result := SASUM(X, incX)
fmt.Println(result)
Output:

3

func SAXPY

func SAXPY(alpha float32, X []float32, incX int, Y []float32, incY int)

Replaces Y by (alpha*X) + Y. Every incX'th and incY'th element is used.

Example
alpha := float32(2)
X := []float32{1, 2, 3}
incX := 1
Y := []float32{4, 5, 6}
incY := 1
SAXPY(alpha, X, incX, Y, incY)
fmt.Println(Y)
Output:

[6 9 12]

func SCASUM

func SCASUM(X []complex64, incX int) float32

Computes the sum of the absolute values of real and imaginary parts of elements in vector X. Every incX'th element is used.

Example
X := []complex64{complex(1, 0), complex(1, 1)}
incX := 1
result := SCASUM(X, incX)
fmt.Println(result)
Output:

3

func SCNRM2

func SCNRM2(X []complex64, incX int) float32

Computes the unitary norm of vector X. Every incX'th element is used.

Example
X := []complex64{complex(1, 0), complex(0, 1)}
incX := 1
result := SCNRM2(X, incX)
fmt.Println(result)
Output:

1.4142135

func SCOPY

func SCOPY(X []float32, incX int, Y []float32, incY int)

Copies X to Y. Every incX'th and incY'th element is used.

Example
X := []float32{1, 666, 2, 666}
incX := 2
Y := []float32{-1, -2}
incY := 1
SCOPY(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[1 666 2 666] [1 2]

func SDOT

func SDOT(X []float32, incX int, Y []float32, incY int) float32

Computes the dot product of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []float32{2, 0}
incX := 1
Y := []float32{3, 5}
incY := 1
result := SDOT(X, incX, Y, incY)
fmt.Println(result)
Output:

6

func SDSDOT

func SDSDOT(alpha float32, X []float32, incX int, Y []float32, incY int) float32

Computes the dot product of vectors X and Y plus an initial value alpha. Every incX'th and incY'th element is used.

Example
alpha := float32(4)
X := []float32{2, 0}
incX := 1
Y := []float32{3, 5}
incY := 1
result := SDSDOT(alpha, X, incX, Y, incY)
fmt.Println(result)
Output:

10

func SGEMM

func SGEMM(transA Transpose, transB Transpose, alpha float32, A [][]float32, B [][]float32, beta float32, C [][]float32)

General matrix-matrix multiplication:

C = alpha*A*B + beta*C

where A and B can optionally be transposed by specifying transA, transB = Trans/NoTrans

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(2, 4)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeFloat32(4, 3)
B[1][0] = 2
B[2][1] = 3

beta := float32(0.0)
C := matrix.MakeFloat32(2, 3)
SGEMM(NoTrans, NoTrans, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[1 -1 0 0] [0 2 0 0]] * [[0 0 0] [2 0 0] [0 3 0] [0 0 0]] = [[-2 0 0] [4 0 0]]

func SGEMV

func SGEMV(transA Transpose, alpha float32, A [][]float32, X []float32, incX int, beta float32, Y []float32, incY int)

Matrix-vector multiplication plus vector with optional matrix transpose. When transA == NoTrans this computes:

Y = alpha*A*X + beta*Y

When transA == Trans this computes:

Y = alpha*(A^T)*X + beta*Y

Every incX'th element of X and incY'th element of Y is used. Matrices must be allocated with MakeFloat32Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeFloat32(2, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []float32{-1, 4, 0}
incX := 1
Y := []float32{0, 0}
incY := 1

alpha := float32(1.0)
beta := float32(0.0)

SGEMV(NoTrans, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)

SGEMV(Trans, alpha, A, Y, incX, beta, X, incY)
fmt.Println(A, "^T*", Y, "=", X)
Output:

[[0 1 2] [3 0 0]] * [-1 4 0] = [4 -3]
[[0 1 2] [3 0 0]] ^T* [4 -3] = [-9 4 8]

func SGER

func SGER(alpha float32, X []float32, incX int, Y []float32, incY int, A [][]float32)

Computes

A = alpha*X*(Y^T) + A

Matrices must be allocated with MakeFloat32Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := float32(1)
X := []float32{1, 2}
incX := 1
Y := []float32{3, 4, 5}
incY := 1
A := matrix.MakeFloat32(2, 3)
SGER(alpha, X, incX, Y, incY, A)
fmt.Println(X, "*", Y, "^T = ", A)
Output:

[1 2] * [3 4 5] ^T =  [[3 4 5] [6 8 10]]

func SNRM2

func SNRM2(X []float32, incX int) float32

Computes the L2 norm (Euclidian length) of vector X. Every incX'th element is used.

Example
X := []float32{1, 666, 1, 666}
incX := 2
result := SNRM2(X, incX)
fmt.Println(result)
Output:

1.4142135

func SROT

func SROT(X []float32, incX int, Y []float32, incY int, c float32, s float32)

Applies a plane rotation to the two vectors X and Y, This routine computes:

[ x_i ]   [ c s ] [ x_i ]
[ y_i ] = [-s c ] [ y_i ]

for all i with strides incX and incY.

Example
X := []float32{1, 0, 1}
incX := 1
Y := []float32{0, 1, 1}
incY := 1
theta := math.Pi / 4
c := float32(math.Cos(theta))
s := float32(math.Sin(theta))
SROT(X, incX, Y, incY, c, s)
fmt.Println(X, Y)
Output:

[0.70710677 0.70710677 1.4142135] [-0.70710677 0.70710677 0]

func SROTG

func SROTG(a float32, b float32) (r, c, s float32)

Constructs a Givens rotation matrix. Computes the values of c and s so that:

[ c  s ] [ a ]    [ r ]
[ -s c ] [ b ] =  [ 0 ].

Returns r, c and s.

Example
a := float32(1)
b := float32(1)
r, c, s := SROTG(a, b)
fmt.Println(r, c, s)
Output:

1.4142135 0.70710677 0.70710677

func SROTM

func SROTM(X []float32, incX int, Y []float32, incY int, P [5]float32)

Applies the modified-Givens rotation of the vectors X and Y:

[ x_i ]   [ h_11 h_12 ] [ x_i ]
[ y_i ] = [ h_21 h_22 ] [ y_i ]

for all i.

H stored in P according to the encoding used by SROTMG.

Example
d1 := float32(1.0)
d2 := float32(1.0)
b1 := float32(1.0)
b2 := float32(1.0)
P := SROTMG(d1, d2, b1, b2)

X := []float32{1}
incX := 1
Y := []float32{1}
incY := 1
SROTM(X, incX, Y, incY, P)
fmt.Println(X, Y)
Output:

[2] [0]

func SROTMG

func SROTMG(d1 float32, d2 float32, b1 float32, b2 float32) [5]float32

Constructs a modified Givens rotation matrix. The input scalars d1, d2, x1 and y1 define a 2-vector [a1 a2]' such that

[ b1 ]   [ d1^{1/2}  0      ] [ x1 ]
[ b2 ] = [   0     d2^{1/2} ] [ y1 ].

This subroutine determines the modified Givens rotation matrix H that transforms y1 and thus a2 to zero. A representation of this matrix is returned:

P[0]: defines the form of matrix H:

-2.0: matrix H contains the identity matrix.
-1.0: matrix H is identical to matrix SH (defined by the remaining values in the vector).
 0.0: H[1,2] and H[2,1] are obtained from matrix SH. The remaining values are both 1.0.
 1.0: H[1,1] and H[2,2] are obtained from matrix SH. H[1,2] is 1.0. H[2,1] is -1.0.

The other elements contain SH:

P[1]	contains SH[1,1].
P[2]	contains SH[2,1].
P[3]	contains SH[1,2].
P[4]	contains SH[2,2].
Example
d1 := float32(1)
d2 := float32(1)
b1 := float32(1)
b2 := float32(1)
result := SROTMG(d1, d2, b1, b2)
fmt.Println(result)
Output:

[1 1 0 0 1]

func SSCAL

func SSCAL(alpha float32, X []float32, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := float32(2)
X := []float32{1, 666, 2, 666}
incX := 2
SSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[2 666 4 666]

func SSWAP

func SSWAP(X []float32, incX int, Y []float32, incY int)

Exchanges the elements of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []float32{1, 666, 2, 666}
incX := 2
Y := []float32{-1, -2}
incY := 1
SSWAP(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[-1 666 -2 666] [1 2]

func SSYMM

func SSYMM(side Side, uplo Uplo, alpha float32, A [][]float32, B [][]float32, beta float32, C [][]float32)

Computes the matrix-matrix product

C = alpha*A*B + beta*C  (for Side==Left)
C = alpha*B*A + beta*C  (for Side==Right)

where the matrix A is symmetric. Only its upper half or lower half is used, specified by uplo=Upper/Lower.

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(3, 3)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeFloat32(3, 2)
B[1][0] = 2
B[2][1] = 3

beta := float32(0.0)
C := matrix.MakeFloat32(3, 2)
SSYMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[1 -1 0] [0 2 0] [0 0 0]] * [[0 0] [2 0] [0 3]] = [[-2 0] [4 0] [0 0]]

func SSYMV

func SSYMV(uplo Uplo, alpha float32, A [][]float32, X []float32, incX int, beta float32, Y []float32, incY int)

Symmetric matrix-vector multiplication:

Y = alpha*A*X + beta*Y

Matrices must be allocated with MakeFloat32Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. A is a symmetric matrix, where uplo (Upper or Lower) determines which part of A is used. Every incX'th element of X and incY'th element of Y is used.

Example
A := matrix.MakeFloat32(2, 2)
A[0][0] = 1
A[1][1] = 2
A[0][1] = 3
A[1][0] = 3

X := []float32{-1, 4}
incX := 1
Y := []float32{0, 0}
incY := 1

alpha := float32(1.0)
beta := float32(0.0)

SSYMV(Upper, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)
Output:

[[1 3] [3 2]] * [-1 4] = [11 5]

func SSYR

func SSYR(uplo Uplo, alpha float32, X []float32, incX int, A [][]float32)

Computes

A = A + alpha*X*(X^T)

A is a symmetric matrix, where uplo (Upper or Lower) determines which part of A is used.

Example
alpha := float32(1)
X := []float32{1, 2}
incX := 1
A := matrix.MakeFloat32(2, 2)
SSYR(Upper, alpha, X, incX, A)
fmt.Println(X, "*", X, "^T = ", A)
Output:

[1 2] * [1 2] ^T =  [[1 2] [0 4]]

func SSYR2

func SSYR2(uplo Uplo, alpha float32, X []float32, incX int, Y []float32, incY int, A [][]float32)

Computes the symmetric rank-2 update

A = \alpha x y^T + \alpha y x^T + A

of the symmetric matrix A. Since the matrix A is symmetric only its upper half or lower half need to be stored. When Uplo is Upper then the upper triangle and diagonal of A are used, and when Uplo is Lower then the lower triangle and diagonal of A are used.

Example
alpha := float32(1)
X := []float32{1, 2}
incX := 1
Y := []float32{2, 4}
incY := 1
A := matrix.MakeFloat32(2, 2)
SSYR2(Upper, alpha, X, incX, Y, incY, A)
fmt.Println(A)
Output:

[[4 8] [0 16]]

func SSYR2K

func SSYR2K(uplo Uplo, trans Transpose, alpha float32, A [][]float32, B [][]float32, beta float32, C [][]float32)

Computes a rank-2k update of the symmetric matrix C:

C = alpha A*(B^T) + alpha B*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*B + alpha (B^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(2, 3)
A[0][0] = 1

B := matrix.MakeFloat32(2, 3)
B[0][0] = 1

beta := float32(0.0)
C := matrix.MakeFloat32(3, 3)

SSYR2K(Upper, Trans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[2 0 0] [0 0 0] [0 0 0]]

func SSYRK

func SSYRK(uplo Uplo, trans Transpose, alpha float32, A [][]float32, beta float32, C [][]float32)

Computes a rank-k update of the symmetric matrix C:

C = alpha A*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(2, 3)
A[0][0] = 1

beta := float32(0.0)
C := matrix.MakeFloat32(3, 3)

SSYRK(Upper, Trans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[1 0 0] [0 0 0] [0 0 0]]

func STRMM

func STRMM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha float32, A [][]float32, B [][]float32)

Computes the matrix-matrix product

B = alpha op(A) B (for Side is Left)
B = alpha B op(A) (for Side is Right)

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(3, 3)
A[1][1] = 1
A[2][2] = 1
B := matrix.MakeFloat32(3, 2)
B[1][1] = 1

STRMM(Left, Upper, NoTrans, NonUnit, alpha, A, B)
fmt.Println(B)
Output:

[[0 0] [0 1] [0 0]]

func STRMV

func STRMV(uplo Uplo, transA Transpose, diag Diag, A [][]float32, X []float32, incX int)

Triangular matrix-vector multiplication plus vector with optional matrix transpose. With uplo == Upper or Lower, the upper or lower triangular part of A is used. diag specifies whether the matrix is unit triangular (). When transA == NoTrans this computes:

X = A*X

When transA == Trans this computes:

X = (A^T)*X

Every incX'th element of X is used. Matrices must be allocated with MakeFloat32Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeFloat32(2, 2)
A[0][0] = 1
A[1][1] = 1

A[0][1] = 2

X := []float32{-1, 4}
incX := 1

fmt.Println(A, "*", X)
STRMV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println("=", X)
Output:

[[1 2] [0 1]] * [-1 4]
= [7 4]

func STRSM

func STRSM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha float32, A [][]float32, B [][]float32)

Computes the inverse-matrix matrix product

B = alpha op(inv(A))B for Side is Left
B = alpha B op(inv(A)) for Side is Right.

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := float32(1.0)
A := matrix.MakeFloat32(3, 3)
B := matrix.MakeFloat32(3, 2)
B[1][1] = 1

STRSM(Left, Upper, NoTrans, Unit, alpha, A, B)
fmt.Println(B)
Output:

[[0 0] [0 1] [0 0]]

func STRSV

func STRSV(uplo Uplo, transA Transpose, diag Diag, A [][]float32, X []float32, incX int)

Computes

inv(op(A)) x for x

where op(A) = A, A^T for TransA = NoTrans, Trans When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of the matrix is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
// solve:
//  x - y = 2
//      y = 3
// for x, y:
A := matrix.MakeFloat32(2, 2)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 1
X := []float32{2, 3}
incX := 1
STRSV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println(X)
Output:

[5 3]

func ZAXPY

func ZAXPY(alpha complex128, X []complex128, incX int, Y []complex128, incY int)

Replaces Y by (alpha*X) + Y. Every incX'th and incY'th element is used.

Example
alpha := complex128(complex(0, 1))
X := []complex128{1, 2, 3}
incX := 1
Y := []complex128{4, 5, 6}
incY := 1
ZAXPY(alpha, X, incX, Y, incY)
fmt.Println(Y)
Output:

[(4+1i) (5+2i) (6+3i)]

func ZCOPY

func ZCOPY(X []complex128, incX int, Y []complex128, incY int)

Copies X to Y. Every incX'th and incY'th element is used.

Example
X := []complex128{1, 666, 2, 666}
incX := 2
Y := []complex128{-1, -2}
incY := 1
ZCOPY(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[(1+0i) (666+0i) (2+0i) (666+0i)] [(1+0i) (2+0i)]

func ZDOTC

func ZDOTC(X []complex128, incX int, Y []complex128, incY int) complex128

Calculates the dot product of the complex conjugate of X with Y. Every incX'th and incY'th element is used.

Example
X := []complex128{complex(1, 0), complex(1, 1)}
incX := 1
Y := []complex128{complex(0, 0), complex(2, 0)}
incY := 1
result := ZDOTC(X, incX, Y, incY)
fmt.Println(result)
Output:

(2-2i)

func ZDOTU

func ZDOTU(X []complex128, incX int, Y []complex128, incY int) complex128

Computes the dot product of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []complex128{complex(1, 0), complex(1, 1)}
incX := 1
Y := []complex128{complex(0, 0), complex(2, 0)}
incY := 1
result := ZDOTU(X, incX, Y, incY)
fmt.Println(result)
Output:

(2+2i)

func ZDSCAL

func ZDSCAL(alpha float64, X []complex128, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := 2.0
X := []complex128{1, 666, 2, 666}
incX := 2
ZDSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[(2+0i) (666+0i) (4+0i) (666+0i)]

func ZGEMM

func ZGEMM(transA Transpose, transB Transpose, alpha complex128, A [][]complex128, B [][]complex128, beta complex128, C [][]complex128)

General matrix-matrix multiplication:

C = alpha*A*B + beta*C

where A and B can optionally be transposed by specifying transA, transB = Trans/NoTrans

Example
alpha := complex128(1.0)
A := matrix.MakeComplex128(2, 4)
A[0][0] = 1
A[0][1] = complex(0, -1)
A[1][1] = 2

B := matrix.MakeComplex128(4, 3)
B[1][0] = 2
B[2][1] = 3

beta := complex128(0.0)
C := matrix.MakeComplex128(2, 3)
ZGEMM(NoTrans, NoTrans, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[(1+0i) (0-1i) (0+0i) (0+0i)] [(0+0i) (2+0i) (0+0i) (0+0i)]] * [[(0+0i) (0+0i) (0+0i)] [(2+0i) (0+0i) (0+0i)] [(0+0i) (3+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]] = [[(0-2i) (0+0i) (0+0i)] [(4+0i) (0+0i) (0+0i)]]

func ZGEMV

func ZGEMV(transA Transpose, alpha complex128, A [][]complex128, X []complex128, incX int, beta complex128, Y []complex128, incY int)

Matrix-vector multiplication plus vector with optional matrix transpose. When transA == NoTrans this computes:

Y = alpha*A*X + beta*Y

When transA == Trans this computes:

Y = alpha*(A^T)*X + beta*Y

Every incX'th element of X and incY'th element of Y is used. Matrices must be allocated with MakeComplex128Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeComplex128(2, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []complex128{-1, 4, 0}
incX := 1
Y := []complex128{0, 0}
incY := 1

alpha := complex(1, 0)
beta := complex(0, 0)

ZGEMV(NoTrans, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A, "*", X, "=", Y)

ZGEMV(Trans, alpha, A, Y, incX, beta, X, incY)
fmt.Println(A, "^T*", Y, "=", X)
Output:

[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)]] * [(-1+0i) (4+0i) (0+0i)] = [(4+0i) (-3+0i)]
[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)]] ^T* [(4+0i) (-3+0i)] = [(-9+0i) (4+0i) (8+0i)]

func ZGERC

func ZGERC(alpha complex128, X []complex128, incX int, Y []complex128, incY int, A [][]complex128)

Computes the conjugate rank-1 update

A = alpha*X*conj(Y^T) + A

Matrices must be allocated with MakeComplex128Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := complex128(1)
X := []complex128{1, 2}
incX := 1
Y := []complex128{complex(0, 3), complex(0, 4), complex(0, 5)}
incY := 1
A := matrix.MakeComplex128(2, 3)
ZGERC(alpha, X, incX, Y, incY, A)
fmt.Println(X, "* conj", Y, "^T = ", A)
Output:

[(1+0i) (2+0i)] * conj [(0+3i) (0+4i) (0+5i)] ^T =  [[(0-3i) (0-4i) (0-5i)] [(0-6i) (0-8i) (0-10i)]]

func ZGERU

func ZGERU(alpha complex128, X []complex128, incX int, Y []complex128, incY int, A [][]complex128)

Computes the rank-1 update:

A = alpha*X*(Y^T) + A

Matrices must be allocated with MakeComplex128Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results. Every incX'th element of X and incY'th element of Y is used.

Example
alpha := complex128(1)
X := []complex128{1, 2}
incX := 1
Y := []complex128{3, 4, 5}
incY := 1
A := matrix.MakeComplex128(2, 3)
ZGERU(alpha, X, incX, Y, incY, A)
fmt.Println(X, "*", Y, "^T = ", A)
Output:

[(1+0i) (2+0i)] * [(3+0i) (4+0i) (5+0i)] ^T =  [[(3+0i) (4+0i) (5+0i)] [(6+0i) (8+0i) (10+0i)]]

func ZHEMM

func ZHEMM(side Side, uplo Uplo, alpha complex128, A [][]complex128, B [][]complex128, beta complex128, C [][]complex128)

Computes the matrix-matrix product and sum

C = alpha A B + beta C (for Side==Left)
C = alpha B A + beta C (for Side==Right)

where the matrix A is hermitian. When Uplo is Upper then the upper triangle and diagonal of A are used, and when Uplo is Lower then the lower triangle and diagonal of A are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := complex128(1.0)
A := matrix.MakeComplex128(2, 2)
A[0][0] = 1
A[0][1] = complex(0, -1)
A[1][1] = 2

B := matrix.MakeComplex128(2, 2)
B[1][0] = 2
B[1][1] = 3

beta := complex128(0.0)
C := matrix.MakeComplex128(2, 2)
ZHEMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(0-2i) (0-3i)] [(4+0i) (6+0i)]]

func ZHEMV

func ZHEMV(uplo Uplo, alpha complex128, A [][]complex128, X []complex128, incX int, beta complex128, Y []complex128, incY int)

Hermitian matrix-vector product:

Y = alpha*X + beta*Y

Matrices must be allocated with MakeComplex128Matrix to ensure contiguous underlying storage, With uplo == Upper or Lower, the upper or lower triangular part of A is used. Every incX'th and incY'th element is used.

Example
A := matrix.MakeComplex128(3, 3)
A[0][1] = 1
A[0][2] = 2
A[1][0] = 3

X := []complex128{-1, 4, 0}
incX := 1
Y := []complex128{0, 0, 0}
incY := 1

alpha := complex128(complex(1, 0))
beta := complex128(complex(0, 0))

ZHEMV(Upper, alpha, A, X, incX, beta, Y, incY)
fmt.Println(A)
Output:

[[(0+0i) (1+0i) (2+0i)] [(3+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func ZHER

func ZHER(uplo Uplo, alpha float64, X []complex128, incX int, A [][]complex128)

Computes the hermitian rank-1 update:

A = alpha*X*conj(X^T) + A

A is a hermitian matrix. Only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := 1.0
X := []complex128{complex(0, 1), 2}
incX := 1
A := matrix.MakeComplex128(2, 2)
ZHER(Upper, alpha, X, incX, A)
fmt.Println(A)
Output:

[[(1+0i) (0+2i)] [(0+0i) (4+0i)]]

func ZHER2

func ZHER2(uplo Uplo, alpha complex128, X []complex128, incX int, Y []complex128, incY int, A [][]complex128)

Computes the hermitian rank-2 update: A = alpha X conj(Y^T) + conj(alpha) Y conj(X^T) + A A is a hermitian matrix. Only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex128(1.0)
X := []complex128{complex(0, 1), 2}
incX := 1
Y := []complex128{2, complex(3, 1)}
incY := 1
A := matrix.MakeComplex128(2, 2)
ZHER2(Upper, alpha, X, incX, Y, incY, A)
fmt.Println(A)
Output:

[[(2+0i) (0+4i)] [(0+0i) (8+0i)]]

func ZHER2K

func ZHER2K(uplo Uplo, trans Transpose, alpha complex128, A [][]complex128, B [][]complex128, beta float64, C [][]complex128)

Computes a rank-2k update of the hermitian matrix C,

C = alpha A B^H + alpha^* B A^H + beta C (for trans==NoTrans)
C = alpha A^H B + alpha^* B^H A + beta C (for Trans==ConjTrans)

Since the matrix C is hermitian only its upper half or lower half need to be stored. When uplo is Upper then the upper triangle and diagonal of C are used, and when uplo is Lower then the lower triangle and diagonal of C are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := complex(2.0, 3.0)
A := matrix.MakeComplex128(2, 2)
A[1][1] = complex(0, 1)
B := matrix.MakeComplex128(2, 2)
B[0][1] = complex(1, 0)
beta := 2.0
C := matrix.MakeComplex128(2, 2)
ZHER2K(Upper, NoTrans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(0+0i) (-3-2i)] [(0+0i) (0+0i)]]

func ZHERK

func ZHERK(uplo Uplo, trans Transpose, alpha float64, A [][]complex128, beta float64, C [][]complex128)

Computes a rank-k update of the hermitian matrix C

C = alpha A A^H + beta C (for trans==NoTrans)
C = alpha A^H A + beta C (for trans==ConjTrans)

Since the matrix C is hermitian only its upper half or lower half need to be stored. When Uplo is Upper then the upper triangle and diagonal of C are used, and when Uplo is Lower then the lower triangle and diagonal of C are used. The imaginary elements of the diagonal are automatically set to zero.

Example
alpha := 1.0
A := matrix.MakeComplex128(2, 2)
A[1][1] = complex(0, 1)
beta := 0.0
C := matrix.MakeComplex128(2, 2)
ZHERK(Upper, ConjTrans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[(0+0i) (0+0i)] [(0+0i) (1+0i)]]

func ZSCAL

func ZSCAL(alpha complex128, X []complex128, incX int)

Multiply X by alpha. Every incX'th element is used.

Example
alpha := complex(0, 1)
X := []complex128{1, 666, 2, 666}
incX := 2
ZSCAL(alpha, X, incX)
fmt.Println(X)
Output:

[(0+1i) (666+0i) (0+2i) (666+0i)]

func ZSWAP

func ZSWAP(X []complex128, incX int, Y []complex128, incY int)

Exchanges the elements of vectors X and Y. Every incX'th and incY'th element is used.

Example
X := []complex128{1, 666, 2, 666}
incX := 2
Y := []complex128{-1, -2}
incY := 1
ZSWAP(X, incX, Y, incY)
fmt.Println(X, Y)
Output:

[(-1+0i) (666+0i) (-2+0i) (666+0i)] [(1+0i) (2+0i)]

func ZSYMM

func ZSYMM(side Side, uplo Uplo, alpha complex128, A [][]complex128, B [][]complex128, beta complex128, C [][]complex128)

Computes the matrix-matrix product

C = alpha*A*B + beta*C  (for Side==Left)
C = alpha*B*A + beta*C  (for Side==Right)

where the matrix A is symmetric. Only its upper half or lower half is used, specified by uplo=Upper/Lower.

Example
alpha := complex128(1)
A := matrix.MakeComplex128(3, 3)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 2

B := matrix.MakeComplex128(3, 2)
B[1][0] = 2
B[2][1] = 3

beta := complex128(0)
C := matrix.MakeComplex128(3, 2)
ZSYMM(Left, Upper, alpha, A, B, beta, C)
fmt.Println(A, "*", B, "=", C)
Output:

[[(1+0i) (-1+0i) (0+0i)] [(0+0i) (2+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]] * [[(0+0i) (0+0i)] [(2+0i) (0+0i)] [(0+0i) (3+0i)]] = [[(-2+0i) (0+0i)] [(4+0i) (0+0i)] [(0+0i) (0+0i)]]

func ZSYR2K

func ZSYR2K(uplo Uplo, trans Transpose, alpha complex128, A [][]complex128, B [][]complex128, beta complex128, C [][]complex128)

Computes a rank-2k update of the symmetric matrix C:

C = alpha A*(B^T) + alpha B*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*B + alpha (B^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex128(1.0)
A := matrix.MakeComplex128(2, 3)
A[0][0] = 1

B := matrix.MakeComplex128(2, 3)
B[0][0] = 1

beta := complex128(0.0)
C := matrix.MakeComplex128(3, 3)

ZSYR2K(Upper, Trans, alpha, A, B, beta, C)
fmt.Println(C)
Output:

[[(2+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func ZSYRK

func ZSYRK(uplo Uplo, trans Transpose, alpha complex128, A [][]complex128, beta complex128, C [][]complex128)

Computes a rank-k update of the symmetric matrix C:

C = alpha A*(A^T) + beta C (for trans==NoTrans)
C = alpha (A^T)*A + beta C (for trans==Trans)

Since the matrix C is symmetric only its upper half or lower half need to be stored, specified by uplo=Upper/Lower.

Example
alpha := complex128(1)
A := matrix.MakeComplex128(2, 3)
A[0][0] = 1

beta := complex128(0)
C := matrix.MakeComplex128(3, 3)

ZSYRK(Upper, Trans, alpha, A, beta, C)
fmt.Println(C)
Output:

[[(1+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)] [(0+0i) (0+0i) (0+0i)]]

func ZTRMM

func ZTRMM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha complex128, A [][]complex128, B [][]complex128)

Computes the matrix-matrix product

B = alpha op(A) B (for Side is Left)
B = alpha B op(A) (for Side is Right)

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := complex128(1)
A := matrix.MakeComplex128(3, 3)
A[1][1] = 1
A[1][2] = complex(0, 1)
B := matrix.MakeComplex128(2, 3)
B[1][2] = 1

ZTRMM(Right, Upper, ConjTrans, NonUnit, alpha, A, B)
fmt.Println(B)
Output:

[[(0+0i) (0+0i) (0+0i)] [(0+0i) (0-1i) (0+0i)]]

func ZTRMV

func ZTRMV(uplo Uplo, transA Transpose, diag Diag, A [][]complex128, X []complex128, incX int)

Triangular matrix-vector multiplication plus vector with optional matrix transpose. With uplo == Upper or Lower, the upper or lower triangular part of A is used. diag specifies whether the matrix is unit triangular (). When transA == NoTrans this computes:

X = A*X

When transA == Trans this computes:

X = (A^T)*X

Every incX'th element of X is used. Matrices must be allocated with MakeComplex128Matrix to ensure contiguous underlying storage, otherwise this function may panic or return unexpected results.

Example
A := matrix.MakeComplex128(2, 2)
A[0][0] = 1
A[1][1] = 1

A[0][1] = 2

X := []complex128{-1, 4}
incX := 1

fmt.Println(A, "*", X)
ZTRMV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println("=", X)
Output:

[[(1+0i) (2+0i)] [(0+0i) (1+0i)]] * [(-1+0i) (4+0i)]
= [(7+0i) (4+0i)]

func ZTRSM

func ZTRSM(side Side, uplo Uplo, transA Transpose, diag Diag, alpha complex128, A [][]complex128, B [][]complex128)

Computes the inverse-matrix matrix product

B = alpha op(inv(A))B for Side is Left
B = alpha B op(inv(A)) for Side is Right.

The matrix A is triangular and op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of A is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
alpha := complex128(1)
A := matrix.MakeComplex128(3, 3)
A[1][1] = 1
A[1][2] = complex(0, 1)
B := matrix.MakeComplex128(2, 3)
B[1][2] = 1

ZTRSM(Right, Upper, ConjTrans, Unit, alpha, A, B)
fmt.Println(B)
Output:

[[(0+0i) (0+0i) (0+0i)] [(0+0i) (0+1i) (1+0i)]]

func ZTRSV

func ZTRSV(uplo Uplo, transA Transpose, diag Diag, A [][]complex128, X []complex128, incX int)

Computes

inv(op(A)) x for x

where op(A) = A, A^T, A^H for TransA = NoTrans, Trans, ConjTrans. When Uplo is Upper then the upper triangle of A is used, and when Uplo is Lower then the lower triangle of A is used. If Diag is NonUnit then the diagonal of the matrix is used, but if Diag is Unit then the diagonal elements of the matrix A are taken as unity and are not referenced.

Example
// solve:
//  x - y = 2
//      y = 3
// for x, y:
A := matrix.MakeComplex128(2, 2)
A[0][0] = 1
A[0][1] = -1
A[1][1] = 1
X := []complex128{2, 3}
incX := 1
ZTRSV(Upper, NoTrans, NonUnit, A, X, incX)
fmt.Println(X)
Output:

[(5+0i) (3+0i)]

Types

type Diag

type Diag uint32 // Used to indicate whether a triangular matrix is unit-diagonal (diagonal elements are all equal to 1).
const (
	NonUnit Diag = C.CblasNonUnit // NonUnit means non unit-diagonal matrix.
	Unit    Diag = C.CblasUnit    //  Unit means unit-diagonal matrix.
)

type Order

type Order uint32 //Indicates whether a matrix is in Row Major or Column Major order.
const (
	RowMajor Order = C.CblasRowMajor // Row major order is the native order for C and Go programs.
	ColMajor Order = C.CblasColMajor // Column major order is native for Fortran.
)

type Side

type Side uint32 // Used to indicate the order of a matrix-matrix multiply.
const (
	Left  Side = C.CblasLeft  // Left means AB matrix-matrix multiply.
	Right Side = C.CblasRight // Right means BA matrix-matrix multiply.
)

type Transpose

type Transpose uint32 // Used to represent transpose operations on a matrix.
const (
	NoTrans   Transpose = C.CblasNoTrans   // NoTrans represents X.
	Trans     Transpose = C.CblasTrans     // Trans represents X^T.
	ConjTrans Transpose = C.CblasConjTrans // ConjTrans represents X^H (=conj(X^T))
)

type Uplo

type Uplo uint32 // Used to indicate which part of a symmetric matrix to use.
const (
	Upper Uplo = C.CblasUpper // Upper means use the upper triangle of the matrix.
	Lower Uplo = C.CblasLower // Lower means use the lower triangle of the matrix.
)

Jump to

Keyboard shortcuts

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