native

package
v0.0.0-...-40ed44c Latest Latest
Warning

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

Go to latest
Published: Feb 23, 2015 License: BSD-3-Clause Imports: 7 Imported by: 0

Documentation

Overview

Package native is a Go implementation of the BLAS API. This implementation panics when the input arguments are invalid as per the standard, for example if a vector increment is zero. 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.

See http://www.netlib.org/lapack/explore-html/d4/de1/_l_i_c_e_n_s_e_source.html for more license information.

Vector arguments are effectively strided slices. They have two input arguments, a number of elements, n, and an increment, incX. The increment specifies the distance between elements of the vector. The actual Go slice may be longer than necessary. The increment may be positive or negative, except in functions with only a single vector argument where the increment may only be positive. If the increment is negative, s[0] is the last element in the slice. Note that this is not the same as counting backward from the end of the slice, as len(s) may be longer than necessary. So, for example, if n = 5 and incX = 3, the elements of s are

[0 * * 1 * * 2 * * 3 * * 4 * * * ...]

where ∗ elements are never accessed. If incX = -3, the same elements are accessed, just in reverse order (4, 3, 2, 1, 0).

Dense matrices are specified by a number of rows, a number of columns, and a stride. The stride specifies the number of entries in the slice between the first element of successive rows. The stride must be at least as large as the number of columns but may be longer.

[a00 ... a0n a0* ... a1stride-1 a21 ... amn am* ... amstride-1]

Thus, dense[i*ld + j] refers to the {i, j}th element of the matrix.

Symmetric and triangular matrices (non-packed) are stored identically to Dense, except that only elements in one triangle of the matrix are accessed.

Packed symmetric and packed triangular matrices are laid out with the entries condensed such that all of the unreferenced elements are removed. So, the upper triangular matrix

[
  1  2  3
  0  4  5
  0  0  6
]

and the lower-triangular matrix

[
  1  0  0
  2  3  0
  4  5  6
]

will both be compacted as [1 2 3 4 5 6]. The (i, j) element of the original dense matrix can be found at element i*n - (i-1)*i/2 + j for upper triangular, and at element i * (i+1) /2 + j for lower triangular.

Banded matrices are laid out in a compact format, constructed by removing the zeros in the rows and aligning the diagonals. For example, the matrix

[
  1  2  3  0  0  0
  4  5  6  7  0  0
  0  8  9 10 11  0
  0  0 12 13 14 15
  0  0  0 16 17 18
  0  0  0  0 19 20
]

implicitly becomes (∗ entries are never accessed)

[
   *  1  2  3
   4  5  6  7
   8  9 10 11
  12 13 14 15
  16 17 18  *
  19 20  *  *
]

which is given to the BLAS routine is [∗ 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) 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

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.

tA and tB specify whether A or B are transposed. A, B, and C are n×n dense matrices.

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 float64, s float64)

Drot applies a plane transformation.

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

func (Implementation) Drotg

func (Implementation) Drotg(a, b float64) (c, s, r, z float64)

Drotg computes the plane rotation

 _    _      _ _       _ _
| c  s |    | a |     | r |
| -s c |  * | b |   = | 0 |
 ‾    ‾      ‾ ‾       ‾ ‾

where

r = ±(a^2 + b^2)
c = a/r, the cosine of the plane rotation
s = b/r, the sine of the plane rotation

NOTE: There is a discrepancy between the refence implementation and the BLAS technical manual regarding the sign for r when a or b are zero. Drotg agrees with the definition in the manual and other common BLAS implementations.

func (Implementation) Drotm

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

Drotm applies the modified Givens rotation to the 2⨉n matrix.

func (Implementation) Drotmg

func (Implementation) Drotmg(d1, d2, x1, y1 float64) (p blas.DrotmParams, rd1, rd2, rx1 float64)

Drotmg computes the modified Givens rotation. See http://www.netlib.org/lapack/explore-html/df/deb/drotmg_8f.html for more details.

func (Implementation) Dsbmv

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

Dsbmv performs

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

func (Implementation) Dspmv(ul blas.Uplo, n int, alpha float64, a []float64, 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, a []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, a []float64)

Dspr2 performs the symmetric rank-2 update

a += alpha * x * y^T + alpha * y * x^T

where a is an n×n symmetric matirx in packed format and x and y are vectors.

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 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, tA 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, tA 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

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 []float64, x []float64, incX int)

Dtpmv computes

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

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, a []float64, x []float64, incX int)

Dtpsv solves

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

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 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 and side == blas.Right

where A is an n×n 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

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

where A is an n×n 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

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

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

Idamax returns the index of the largest element of x. If there are multiple such indices the earliest is returned. Idamax returns -1 if incX is negative or if n == 0.

Jump to

Keyboard shortcuts

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