gonum: gonum.org/v1/gonum/mat Index | Examples | Files

package mat

import "gonum.org/v1/gonum/mat"

Package mat provides implementations of float64 and complex128 matrix structures and linear algebra operations on them.

Overview

This section provides a quick overview of the mat package. The following sections provide more in depth commentary.

mat provides:

- Interfaces for Matrix classes (Matrix, Symmetric, Triangular)
- Concrete implementations (Dense, SymDense, TriDense)
- Methods and functions for using matrix data (Add, Trace, SymRankOne)
- Types for constructing and using matrix factorizations (QR, LU)
- The complementary types for complex matrices, CMatrix, CSymDense, etc.

A matrix may be constructed through the corresponding New function. If no backing array is provided the matrix will be initialized to all zeros.

// Allocate a zeroed real matrix of size 3×5
zero := mat.NewDense(3, 5, nil)

If a backing data slice is provided, the matrix will have those elements. Matrices are all stored in row-major format.

// Generate a 6×6 matrix of random values.
data := make([]float64, 36)
for i := range data {
	data[i] = rand.NormFloat64()
}
a := mat.NewDense(6, 6, data)

Operations involving matrix data are implemented as functions when the values of the matrix remain unchanged

tr := mat.Trace(a)

and are implemented as methods when the operation modifies the receiver.

zero.Copy(a)

Receivers must be the correct size for the matrix operations, otherwise the operation will panic. As a special case for convenience, a zero-value matrix will be modified to have the correct size, allocating data if necessary.

var c mat.Dense // construct a new zero-sized matrix
c.Mul(a, a)     // c is automatically adjusted to be 6×6

Zero-value of a matrix

A zero-value matrix is either the Go language definition of a zero-value or is a zero-sized matrix with zero-length stride. Matrix implementations may have a Reset method to revert the receiver into a zero-valued matrix and an IsZero method that returns whether the matrix is zero-valued. So the following will all result in a zero-value matrix.

- var a mat.Dense
- a := NewDense(0, 0, make([]float64, 0, 100))
- a.Reset()

A zero-value matrix can not be sliced even if it does have an adequately sized backing data slice, but can be expanded using its Grow method if it exists.

The Matrix Interfaces

The Matrix interface is the common link between the concrete types of real matrices, The Matrix interface is defined by three functions: Dims, which returns the dimensions of the Matrix, At, which returns the element in the specified location, and T for returning a Transpose (discussed later). All of the concrete types can perform these behaviors and so implement the interface. Methods and functions are designed to use this interface, so in particular the method

func (m *Dense) Mul(a, b Matrix)

constructs a *Dense from the result of a multiplication with any Matrix types, not just *Dense. Where more restrictive requirements must be met, there are also the Symmetric and Triangular interfaces. For example, in

func (s *SymDense) AddSym(a, b Symmetric)

the Symmetric interface guarantees a symmetric result.

The CMatrix interface plays the same role for complex matrices. The difference is that the CMatrix type has the H method instead T, for returning the conjugate transpose.

(Conjugate) Transposes

The T method is used for transposition on real matrices, and H is used for conjugate transposition on complex matrices. For example, c.Mul(a.T(), b) computes c = a^T * b. The mat types implement this method implicitly — see the Transpose and Conjugate types for more details. Note that some operations have a transpose as part of their definition, as in *SymDense.SymOuterK.

Matrix Factorization

Matrix factorizations, such as the LU decomposition, typically have their own specific data storage, and so are each implemented as a specific type. The factorization can be computed through a call to Factorize

var lu mat.LU
lu.Factorize(a)

The elements of the factorization can be extracted through methods on the factorized type, i.e. *LU.UTo. The factorization types can also be used directly, as in *Dense.SolveCholesky. Some factorizations can be updated directly, without needing to update the original matrix and refactorize, as in *LU.RankOne.

BLAS and LAPACK

BLAS and LAPACK are the standard APIs for linear algebra routines. Many operations in mat are implemented using calls to the wrapper functions in gonum/blas/blas64 and gonum/lapack/lapack64 and their complex equivalents. By default, blas64 and lapack64 call the native Go implementations of the routines. Alternatively, it is possible to use C-based implementations of the APIs through the respective cgo packages and "Use" functions. The Go implementation of LAPACK (used by default) makes calls through blas64, so if a cgo BLAS implementation is registered, the lapack64 calls will be partially executed in Go and partially executed in C.

Type Switching

The Matrix abstraction enables efficiency as well as interoperability. Go's type reflection capabilities are used to choose the most efficient routine given the specific concrete types. For example, in

c.Mul(a, b)

if a and b both implement RawMatrixer, that is, they can be represented as a blas64.General, blas64.Gemm (general matrix multiplication) is called, while instead if b is a RawSymmetricer blas64.Symm is used (general-symmetric multiplication), and if b is a *VecDense blas64.Gemv is used.

There are many possible type combinations and special cases. No specific guarantees are made about the performance of any method, and in particular, note that an abstract matrix type may be copied into a concrete type of the corresponding value. If there are specific special cases that are needed, please submit a pull-request or file an issue.

Invariants

Matrix input arguments to functions are never directly modified. If an operation changes Matrix data, the mutated matrix will be the receiver of a function.

For convenience, a matrix may be used as both a receiver and as an input, e.g.

a.Pow(a, 6)
v.SolveVec(a.T(), v)

though in many cases this will cause an allocation (see Element Aliasing). An exception to this rule is Copy, which does not allow a.Copy(a.T()).

Element Aliasing

Most methods in mat modify receiver data. It is forbidden for the modified data region of the receiver to overlap the used data area of the input arguments. The exception to this rule is when the method receiver is equal to one of the input arguments, as in the a.Pow(a, 6) call above, or its implicit transpose.

This prohibition is to help avoid subtle mistakes when the method needs to read from and write to the same data region. There are ways to make mistakes using the mat API, and mat functions will detect and complain about those. There are many ways to make mistakes by excursion from the mat API via interaction with raw matrix values.

If you need to read the rest of this section to understand the behavior of your program, you are being clever. Don't be clever. If you must be clever, blas64 and lapack64 may be used to call the behavior directly.

mat will use the following rules to detect overlap between the receiver and one of the inputs:

- the input implements one of the Raw methods, and
- the address ranges of the backing data slices overlap, and
- the strides differ or there is an overlap in the used data elements.

If such an overlap is detected, the method will panic.

The following cases will not panic:

- the data slices do not overlap,
- there is pointer identity between the receiver and input values after
  the value has been untransposed if necessary.

mat will not attempt to detect element overlap if the input does not implement a Raw method. Method behavior is undefined if there is undetected overlap.

Index

Examples

Package Files

band.go cholesky.go cmatrix.go consts.go dense.go dense_arithmetic.go doc.go eigen.go errors.go format.go gsvd.go hogsvd.go index_no_bound_checks.go inner.go io.go lq.go lu.go matrix.go offset.go pool.go product.go qr.go shadow.go solve.go svd.go symband.go symmetric.go triangular.go vector.go

Constants

const (
    // CondNorm is the matrix norm used for computing the condition number by routines
    // in the matrix packages.
    CondNorm = lapack.MaxRowSum

    // CondNormTrans is the norm used to compute on A^T to get the same result as
    // computing CondNorm on A.
    CondNormTrans = lapack.MaxColumnSum
)
const ConditionTolerance = 1e16

ConditionTolerance is the tolerance limit of the condition number. If the condition number is above this value, the matrix is considered singular.

Variables

var (
    ErrIndexOutOfRange     = Error{"matrix: index out of range"}
    ErrRowAccess           = Error{"matrix: row index out of range"}
    ErrColAccess           = Error{"matrix: column index out of range"}
    ErrVectorAccess        = Error{"matrix: vector index out of range"}
    ErrZeroLength          = Error{"matrix: zero length in matrix definition"}
    ErrRowLength           = Error{"matrix: row length mismatch"}
    ErrColLength           = Error{"matrix: col length mismatch"}
    ErrSquare              = Error{"matrix: expect square matrix"}
    ErrNormOrder           = Error{"matrix: invalid norm order for matrix"}
    ErrSingular            = Error{"matrix: matrix is singular"}
    ErrShape               = Error{"matrix: dimension mismatch"}
    ErrIllegalStride       = Error{"matrix: illegal stride"}
    ErrPivot               = Error{"matrix: malformed pivot list"}
    ErrTriangle            = Error{"matrix: triangular storage mismatch"}
    ErrTriangleSet         = Error{"matrix: triangular set out of bounds"}
    ErrBandSet             = Error{"matrix: band set out of bounds"}
    ErrSliceLengthMismatch = Error{"matrix: input slice length mismatch"}
    ErrNotPSD              = Error{"matrix: input not positive symmetric definite"}
    ErrFailedEigen         = Error{"matrix: eigendecomposition not successful"}
)

func Col Uses

func Col(dst []float64, j int, a Matrix) []float64

Col copies the elements in the jth column of the matrix into the slice dst. The length of the provided slice must equal the number of rows, unless the slice is nil in which case a new slice is first allocated.

func Cond Uses

func Cond(a Matrix, norm float64) float64

Cond returns the condition number of the given matrix under the given norm. The condition number must be based on the 1-norm, 2-norm or ∞-norm. Cond will panic with matrix.ErrShape if the matrix has zero size.

BUG(btracey): The computation of the 1-norm and ∞-norm for non-square matrices is innacurate, although is typically the right order of magnitude. See https://github.com/xianyi/OpenBLAS/issues/636. While the value returned will change with the resolution of this bug, the result from Cond will match the condition number used internally.

func Det Uses

func Det(a Matrix) float64

Det returns the determinant of the matrix a. In many expressions using LogDet will be more numerically stable.

func Dot Uses

func Dot(a, b Vector) float64

Dot returns the sum of the element-wise product of a and b. Dot panics if the matrix sizes are unequal.

func Equal Uses

func Equal(a, b Matrix) bool

Equal returns whether the matrices a and b have the same size and are element-wise equal.

func EqualApprox Uses

func EqualApprox(a, b Matrix, epsilon float64) bool

EqualApprox returns whether the matrices a and b have the same size and contain all equal elements with tolerance for element-wise equality specified by epsilon. Matrices with non-equal shapes are not equal.

func Formatted Uses

func Formatted(m Matrix, options ...FormatOption) fmt.Formatter

Formatted returns a fmt.Formatter for the matrix m using the given options.

Code:

a := mat.NewDense(3, 3, []float64{1, 2, 3, 0, 4, 5, 0, 0, 6})

// Create a matrix formatting value with a prefix and calculating each column
// width individually...
fa := mat.Formatted(a, mat.Prefix("    "), mat.Squeeze())

// and then print with and without zero value elements.
fmt.Printf("with all values:\na = %v\n\n", fa)
fmt.Printf("with only non-zero values:\na = % v\n\n", fa)

// Modify the matrix...
a.Set(0, 2, 0)

// and print it without zero value elements.
fmt.Printf("after modification with only non-zero values:\na = % v\n\n", fa)

// Modify the matrix again...
a.Set(0, 2, 123.456)

// and print it using scientific notation for large exponents.
fmt.Printf("after modification with scientific notation:\na = %.2g\n\n", fa)
// See golang.org/pkg/fmt/ floating-point verbs for a comprehensive list.

Output:

with all values:
a = ⎡1  2  3⎤
    ⎢0  4  5⎥
    ⎣0  0  6⎦

with only non-zero values:
a = ⎡1  2  3⎤
    ⎢.  4  5⎥
    ⎣.  .  6⎦

after modification with only non-zero values:
a = ⎡1  2  .⎤
    ⎢.  4  5⎥
    ⎣.  .  6⎦

after modification with scientific notation:
a = ⎡1  2  1.2e+02⎤
    ⎢0  4        5⎥
    ⎣0  0        6⎦

func Inner Uses

func Inner(x *VecDense, A Matrix, y *VecDense) float64

Inner computes the generalized inner product

x^T A y

between vectors x and y with matrix A. This is only a true inner product if A is symmetric positive definite, though the operation works for any matrix A.

Inner panics if x.Len != m or y.Len != n when A is an m x n matrix.

func LogDet Uses

func LogDet(a Matrix) (det float64, sign float64)

LogDet returns the log of the determinant and the sign of the determinant for the matrix that has been factorized. Numerical stability in product and division expressions is generally improved by working in log space.

func Max Uses

func Max(a Matrix) float64

Max returns the largest element value of the matrix A. Max will panic with matrix.ErrShape if the matrix has zero size.

func Maybe Uses

func Maybe(fn func()) (err error)

Maybe will recover a panic with a type mat.Error from fn, and return this error as the Err field of an ErrorStack. The stack trace for the panicking function will be recovered and placed in the StackTrace field. Any other error is re-panicked.

func MaybeComplex Uses

func MaybeComplex(fn func() complex128) (f complex128, err error)

MaybeComplex will recover a panic with a type mat.Error from fn, and return this error as the Err field of an ErrorStack. The stack trace for the panicking function will be recovered and placed in the StackTrace field. Any other error is re-panicked.

func MaybeFloat Uses

func MaybeFloat(fn func() float64) (f float64, err error)

MaybeFloat will recover a panic with a type mat.Error from fn, and return this error as the Err field of an ErrorStack. The stack trace for the panicking function will be recovered and placed in the StackTrace field. Any other error is re-panicked.

func Min Uses

func Min(a Matrix) float64

Min returns the smallest element value of the matrix A. Min will panic with matrix.ErrShape if the matrix has zero size.

func Norm Uses

func Norm(a Matrix, norm float64) float64

Norm returns the specified (induced) norm of the matrix a. See https://en.wikipedia.org/wiki/Matrix_norm for the definition of an induced norm.

Valid norms are:

  1 - The maximum absolute column sum
  2 - Frobenius norm, the square root of the sum of the squares of the elements.
Inf - The maximum absolute row sum.

Norm will panic with ErrNormOrder if an illegal norm order is specified and with matrix.ErrShape if the matrix has zero size.

func Row Uses

func Row(dst []float64, i int, a Matrix) []float64

Row copies the elements in the ith row of the matrix into the slice dst. The length of the provided slice must equal the number of columns, unless the slice is nil in which case a new slice is first allocated.

func Sum Uses

func Sum(a Matrix) float64

Sum returns the sum of the elements of the matrix.

func Trace Uses

func Trace(a Matrix) float64

Trace returns the trace of the matrix. Trace will panic if the matrix is not square.

type BandDense Uses

type BandDense struct {
    // contains filtered or unexported fields
}

BandDense represents a band matrix in dense storage format.

func NewBandDense Uses

func NewBandDense(r, c, kl, ku int, data []float64) *BandDense

NewBandDense creates a new Band matrix with r rows and c columns. If data == nil, a new slice is allocated for the backing slice. If len(data) == min(r, c+kl)*(kl+ku+1), data is used as the backing slice, and changes to the elements of the returned BandDense will be reflected in data. If neither of these is true, NewBandDense will panic. kl must be at least zero and less r, and ku must be at least zero and less than c, otherwise NewBandDense will panic.

The data must be arranged in row-major order constructed by removing the zeros from the rows outside the band 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

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 passed to NewBandDense as []float64{*, 1, 2, 3, 4, ...} with kl=1 and ku=2. Only the values in the band portion of the matrix are used.

func NewDiagonalRect Uses

func NewDiagonalRect(r, c int, data []float64) *BandDense

NewDiagonalRect is a convenience function that returns a diagonal matrix represented by a BandDense. The length of data must be min(r, c) otherwise NewDiagonalRect will panic.

func (*BandDense) At Uses

func (b *BandDense) At(i, j int) float64

At returns the element at row i, column j.

func (*BandDense) Bandwidth Uses

func (b *BandDense) Bandwidth() (kl, ku int)

Bandwidth returns the upper and lower bandwidths of the matrix.

func (*BandDense) Dims Uses

func (b *BandDense) Dims() (r, c int)

Dims returns the number of rows and columns in the matrix.

func (*BandDense) DoColNonZero Uses

func (b *BandDense) DoColNonZero(j int, fn func(i, j int, v float64))

DoColNonZero calls the function fn for each of the non-zero elements of column j of b. The function fn takes a row/column index and the element value of b at (i, j).

func (*BandDense) DoNonZero Uses

func (b *BandDense) DoNonZero(fn func(i, j int, v float64))

DoNonZero calls the function fn for each of the non-zero elements of b. The function fn takes a row/column index and the element value of b at (i, j).

func (*BandDense) DoRowNonZero Uses

func (b *BandDense) DoRowNonZero(i int, fn func(i, j int, v float64))

DoRowNonZero calls the function fn for each of the non-zero elements of row i of b. The function fn takes a row/column index and the element value of b at (i, j).

func (*BandDense) RawBand Uses

func (b *BandDense) RawBand() blas64.Band

RawBand returns the underlying blas64.Band used by the receiver. Changes to elements in the receiver following the call will be reflected in returned blas64.Band.

func (*BandDense) SetBand Uses

func (b *BandDense) SetBand(i, j int, v float64)

SetBand sets the element at row i, column j to the value v. It panics if the location is outside the appropriate region of the matrix.

func (*BandDense) T Uses

func (b *BandDense) T() Matrix

T performs an implicit transpose by returning the receiver inside a Transpose.

func (*BandDense) TBand Uses

func (b *BandDense) TBand() Banded

TBand performs an implicit transpose by returning the receiver inside a TransposeBand.

type BandWidther Uses

type BandWidther interface {
    BandWidth() (k1, k2 int)
}

A BandWidther represents a banded matrix and can return the left and right half-bandwidths, k1 and k2.

type Banded Uses

type Banded interface {
    Matrix
    // Bandwidth returns the lower and upper bandwidth values for
    // the matrix. The total bandwidth of the matrix is kl+ku+1.
    Bandwidth() (kl, ku int)

    // TBand is the equivalent of the T() method in the Matrix
    // interface but guarantees the transpose is of banded type.
    TBand() Banded
}

Banded is a band matrix representation.

type CMatrix Uses

type CMatrix interface {
    // Dims returns the dimensions of a Matrix.
    Dims() (r, c int)

    // At returns the value of a matrix element at row i, column j.
    // It will panic if i or j are out of bounds for the matrix.
    At(i, j int) complex128

    // H returns the conjugate transpose of the Matrix. Whether H
    // returns a copy of the underlying data is implementation dependent.
    // This method may be implemented using the Conjugate type, which
    // provides an implicit matrix conjugate transpose.
    H() CMatrix
}

CMatrix is the basic matrix interface type for complex matrices.

type Cholesky Uses

type Cholesky struct {
    // contains filtered or unexported fields
}

Cholesky is a type for creating and using the Cholesky factorization of a symmetric positive definite matrix.

Cholesky methods may only be called on a value that has been successfully initialized by a call to Factorize that has returned true. Calls to methods of an unsuccessful Cholesky factorization will panic.

Code:

// Construct a symmetric positive definite matrix.
tmp := mat.NewDense(4, 4, []float64{
    2, 6, 8, -4,
    1, 8, 7, -2,
    2, 2, 1, 7,
    8, -2, -2, 1,
})
var a mat.SymDense
a.SymOuterK(1, tmp)

fmt.Printf("a = %0.4v\n", mat.Formatted(&a, mat.Prefix("    ")))

// Compute the cholesky factorization.
var chol mat.Cholesky
if ok := chol.Factorize(&a); !ok {
    fmt.Println("a matrix is not positive semi-definite.")
}

// Find the determinant.
fmt.Printf("\nThe determinant of a is %0.4g\n\n", chol.Det())

// Use the factorization to solve the system of equations a * x = b.
b := mat.NewVecDense(4, []float64{1, 2, 3, 4})
var x mat.VecDense
if err := chol.SolveVec(&x, b); err != nil {
    fmt.Println("Matrix is near singular: ", err)
}
fmt.Println("Solve a * x = b")
fmt.Printf("x = %0.4v\n", mat.Formatted(&x, mat.Prefix("    ")))

// Extract the factorization and check that it equals the original matrix.
t := chol.LTo(nil)
var test mat.Dense
test.Mul(t, t.T())
fmt.Println()
fmt.Printf("L * L^T = %0.4v\n", mat.Formatted(&a, mat.Prefix("          ")))

Output:

a = ⎡120  114   -4  -16⎤
    ⎢114  118   11  -24⎥
    ⎢ -4   11   58   17⎥
    ⎣-16  -24   17   73⎦

The determinant of a is 1.543e+06

Solve a * x = b
x = ⎡  -0.239⎤
    ⎢  0.2732⎥
    ⎢-0.04681⎥
    ⎣  0.1031⎦

L * L^T = ⎡120  114   -4  -16⎤
          ⎢114  118   11  -24⎥
          ⎢ -4   11   58   17⎥
          ⎣-16  -24   17   73⎦

func (*Cholesky) Clone Uses

func (c *Cholesky) Clone(chol *Cholesky)

Clone makes a copy of the input Cholesky into the receiver, overwriting the previous value of the receiver. Clone does not place any restrictions on receiver shape. Clone panics if the input Cholesky is not the result of a valid decomposition.

func (*Cholesky) Cond Uses

func (c *Cholesky) Cond() float64

Cond returns the condition number of the factorized matrix.

func (*Cholesky) Det Uses

func (c *Cholesky) Det() float64

Det returns the determinant of the matrix that has been factorized.

func (*Cholesky) ExtendVecSym Uses

func (chol *Cholesky) ExtendVecSym(a *Cholesky, v Vector) (ok bool)

ExtendVecSym computes the Cholesky decomposition of the original matrix A, whose Cholesky decomposition is in a, extended by a the n×1 vector v according to

[A  w]
[w' k]

where k = v[n-1] and w = v[:n-1]. The result is stored into the receiver. In order for the updated matrix to be positive definite, it must be the case that k > w' A^-1 w. If this condition does not hold then ExtendVecSym will return false and the receiver will not be updated.

ExtendVecSym will panic if v.Len() != a.Size()+1 or if a does not contain a valid decomposition.

func (*Cholesky) Factorize Uses

func (c *Cholesky) Factorize(a Symmetric) (ok bool)

Factorize calculates the Cholesky decomposition of the matrix A and returns whether the matrix is positive definite. If Factorize returns false, the factorization must not be used.

func (*Cholesky) InverseTo Uses

func (c *Cholesky) InverseTo(s *SymDense) error

InverseTo computes the inverse of the matrix represented by its Cholesky factorization and stores the result into s. If the factorized matrix is ill-conditioned, a Condition error will be returned. Note that matrix inversion is numerically unstable, and should generally be avoided where possible, for example by using the Solve routines.

func (*Cholesky) LTo Uses

func (c *Cholesky) LTo(dst *TriDense) *TriDense

LTo extracts the n×n lower triangular matrix L from a Cholesky decomposition into dst and returns the result. If dst is nil a new TriDense is allocated.

A = L * L^T.

func (*Cholesky) LogDet Uses

func (c *Cholesky) LogDet() float64

LogDet returns the log of the determinant of the matrix that has been factorized.

func (*Cholesky) RawU Uses

func (c *Cholesky) RawU() Triangular

RawU returns the Triangular matrix used to store the Cholesky decomposition of the original matrix A. The returned matrix should not be modified. If it is modified, the decomposition is invalid and should not be used.

func (*Cholesky) Reset Uses

func (c *Cholesky) Reset()

Reset resets the factorization so that it can be reused as the receiver of a dimensionally restricted operation.

func (*Cholesky) Scale Uses

func (c *Cholesky) Scale(f float64, orig *Cholesky)

Scale multiplies the original matrix A by a positive constant using its Cholesky decomposition, storing the result in-place into the receiver. That is, if the original Cholesky factorization is

U^T * U = A

the updated factorization is

U'^T * U' = f A = A'

Scale panics if the constant is non-positive, or if the receiver is non-zero and is of a different Size from the input.

func (*Cholesky) SetFromU Uses

func (c *Cholesky) SetFromU(t *TriDense)

SetFromU sets the Cholesky decomposition from the given triangular matrix. SetFromU panics if t is not upper triangular. Note that t is copied into, not stored inside, the receiver.

func (*Cholesky) Size Uses

func (c *Cholesky) Size() int

Size returns the dimension of the factorized matrix.

func (*Cholesky) Solve Uses

func (c *Cholesky) Solve(m *Dense, b Matrix) error

Solve finds the matrix m that solves A * m = b where A is represented by the Cholesky decomposition, placing the result in m.

func (*Cholesky) SolveChol Uses

func (a *Cholesky) SolveChol(m *Dense, b *Cholesky) error

SolveChol finds the matrix m that solves A * m = B where A and B are represented by their Cholesky decompositions a and b, placing the result in the receiver.

func (*Cholesky) SolveVec Uses

func (c *Cholesky) SolveVec(v, b *VecDense) error

SolveVec finds the vector v that solves A * v = b where A is represented by the Cholesky decomposition, placing the result in v.

func (*Cholesky) SymRankOne Uses

func (c *Cholesky) SymRankOne(orig *Cholesky, alpha float64, x *VecDense) (ok bool)

SymRankOne performs a rank-1 update of the original matrix A and refactorizes its Cholesky factorization, storing the result into the receiver. That is, if in the original Cholesky factorization

U^T * U = A,

in the updated factorization

U'^T * U' = A + alpha * x * x^T = A'.

Note that when alpha is negative, the updating problem may be ill-conditioned and the results may be inaccurate, or the updated matrix A' may not be positive definite and not have a Cholesky factorization. SymRankOne returns whether the updated matrix A' is positive definite.

SymRankOne updates a Cholesky factorization in O(n²) time. The Cholesky factorization computation from scratch is O(n³).

Code:

a := mat.NewSymDense(4, []float64{
    1, 1, 1, 1,
    0, 2, 3, 4,
    0, 0, 6, 10,
    0, 0, 0, 20,
})
fmt.Printf("A = %0.4v\n", mat.Formatted(a, mat.Prefix("    ")))

// Compute the Cholesky factorization.
var chol mat.Cholesky
if ok := chol.Factorize(a); !ok {
    fmt.Println("matrix a is not positive definite.")
}

x := mat.NewVecDense(4, []float64{0, 0, 0, 1})
fmt.Printf("\nx = %0.4v\n", mat.Formatted(x, mat.Prefix("    ")))

// Rank-1 update the factorization.
chol.SymRankOne(&chol, 1, x)
// Rank-1 update the matrix a.
a.SymRankOne(a, 1, x)

au := chol.ToSym(nil)

// Print the matrix that was updated directly.
fmt.Printf("\nA' =        %0.4v\n", mat.Formatted(a, mat.Prefix("            ")))
// Print the matrix recovered from the factorization.
fmt.Printf("\nU'^T * U' = %0.4v\n", mat.Formatted(au, mat.Prefix("            ")))

Output:

A = ⎡ 1   1   1   1⎤
    ⎢ 1   2   3   4⎥
    ⎢ 1   3   6  10⎥
    ⎣ 1   4  10  20⎦

x = ⎡0⎤
    ⎢0⎥
    ⎢0⎥
    ⎣1⎦

A' =        ⎡ 1   1   1   1⎤
            ⎢ 1   2   3   4⎥
            ⎢ 1   3   6  10⎥
            ⎣ 1   4  10  21⎦

U'^T * U' = ⎡ 1   1   1   1⎤
            ⎢ 1   2   3   4⎥
            ⎢ 1   3   6  10⎥
            ⎣ 1   4  10  21⎦

func (*Cholesky) ToSym Uses

func (c *Cholesky) ToSym(dst *SymDense) *SymDense

ToSym reconstructs the original positive definite matrix given its Cholesky decomposition into dst and returns the result. If dst is nil a new SymDense is allocated.

func (*Cholesky) UTo Uses

func (c *Cholesky) UTo(dst *TriDense) *TriDense

UTo extracts the n×n upper triangular matrix U from a Cholesky decomposition into dst and returns the result. If dst is nil a new TriDense is allocated.

A = U^T * U.

type Cloner Uses

type Cloner interface {
    Clone(a Matrix)
}

A Cloner can make a copy of a into the receiver, overwriting the previous value of the receiver. The clone operation does not make any restriction on shape and will not cause shadowing.

type ColNonZeroDoer Uses

type ColNonZeroDoer interface {
    DoColNonZero(j int, fn func(i, j int, v float64))
}

A ColNonZeroDoer can call a function for each non-zero element of a column of the receiver. The parameters of the function are the element indices and its value.

type ColViewer Uses

type ColViewer interface {
    ColView(j int) Vector
}

A ColViewer can return a Vector reflecting a column that is backed by the matrix data. The Vector returned will have length equal to the number of rows.

type Condition Uses

type Condition float64

Condition is the condition number of a matrix. The condition number is defined as |A| * |A^-1|.

One important use of Condition is during linear solve routines (finding x such that A * x = b). The condition number of A indicates the accuracy of the computed solution. A Condition error will be returned if the condition number of A is sufficiently large. If A is exactly singular to working precision, Condition == ∞, and the solve algorithm may have completed early. If Condition is large and finite the solve algorithm will be performed, but the computed solution may be innacurate. Due to the nature of finite precision arithmetic, the value of Condition is only an approximate test of singularity.

func (Condition) Error Uses

func (c Condition) Error() string

type Conjugate Uses

type Conjugate struct {
    CMatrix CMatrix
}

Conjugate is a type for performing an implicit matrix conjugate transpose. It implements the Matrix interface, returning values from the conjugate transpose of the matrix within.

func (Conjugate) At Uses

func (t Conjugate) At(i, j int) complex128

At returns the value of the element at row i and column j of the transposed matrix, that is, row j and column i of the Matrix field.

func (Conjugate) Dims Uses

func (t Conjugate) Dims() (r, c int)

Dims returns the dimensions of the transposed matrix. The number of rows returned is the number of columns in the Matrix field, and the number of columns is the number of rows in the Matrix field.

func (Conjugate) H Uses

func (t Conjugate) H() CMatrix

H performs an implicit conjugate transpose by returning the Matrix field.

func (Conjugate) Unconjugate Uses

func (t Conjugate) Unconjugate() CMatrix

Unconjugate returns the Matrix field.

type Copier Uses

type Copier interface {
    Copy(a Matrix) (r, c int)
}

A Copier can make a copy of elements of a into the receiver. The submatrix copied starts at row and column 0 and has dimensions equal to the minimum dimensions of the two matrices. The number of row and columns copied is returned. Copy will copy from a source that aliases the receiver unless the source is transposed; an aliasing transpose copy will panic with the exception for a special case when the source data has a unitary increment or stride.

type Dense Uses

type Dense struct {
    // contains filtered or unexported fields
}

Dense is a dense matrix representation.

func DenseCopyOf Uses

func DenseCopyOf(a Matrix) *Dense

DenseCopyOf returns a newly allocated copy of the elements of a.

func NewDense Uses

func NewDense(r, c int, data []float64) *Dense

NewDense creates a new Dense matrix with r rows and c columns. If data == nil, a new slice is allocated for the backing slice. If len(data) == r*c, data is used as the backing slice, and changes to the elements of the returned Dense will be reflected in data. If neither of these is true, NewDense will panic.

The data must be arranged in row-major order, i.e. the (i*c + j)-th element in the data slice is the {i, j}-th element in the matrix.

func (*Dense) Add Uses

func (m *Dense) Add(a, b Matrix)

Add adds a and b element-wise, placing the result in the receiver. Add will panic if the two matrices do not have the same shape.

Code:

// Initialize two matrices, a and b.
a := mat.NewDense(2, 2, []float64{
    1, 0,
    1, 0,
})
b := mat.NewDense(2, 2, []float64{
    0, 1,
    0, 1,
})

// Add a and b, placing the result into c.
// Notice that the size is automatically adjusted
// when the receiver has zero size.
var c mat.Dense
c.Add(a, b)

// Print the result using the formatter.
fc := mat.Formatted(&c, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("c = %v", fc)

Output:

c = ⎡1  1⎤
    ⎣1  1⎦

func (*Dense) Apply Uses

func (m *Dense) Apply(fn func(i, j int, v float64) float64, a Matrix)

Apply applies the function fn to each of the elements of a, placing the resulting matrix in the receiver. The function fn takes a row/column index and element value and returns some function of that tuple.

func (*Dense) At Uses

func (m *Dense) At(i, j int) float64

At returns the element at row i, column j.

func (*Dense) Augment Uses

func (m *Dense) Augment(a, b Matrix)

Augment creates the augmented matrix of a and b, where b is placed in the greater indexed columns. Augment will panic if the two input matrices do not have the same number of rows or the constructed augmented matrix is not the same shape as the receiver.

func (*Dense) Caps Uses

func (m *Dense) Caps() (r, c int)

Caps returns the number of rows and columns in the backing matrix.

func (*Dense) Clone Uses

func (m *Dense) Clone(a Matrix)

Clone makes a copy of a into the receiver, overwriting the previous value of the receiver. The clone operation does not make any restriction on shape and will not cause shadowing.

See the Cloner interface for more information.

func (*Dense) ColView Uses

func (m *Dense) ColView(j int) Vector

ColView returns a Vector reflecting the column j, backed by the matrix data.

See ColViewer for more information.

func (*Dense) Copy Uses

func (m *Dense) Copy(a Matrix) (r, c int)

Copy makes a copy of elements of a into the receiver. It is similar to the built-in copy; it copies as much as the overlap between the two matrices and returns the number of rows and columns it copied. If a aliases the receiver and is a transposed Dense or VecDense, with a non-unitary increment, Copy will panic.

See the Copier interface for more information.

func (*Dense) Dims Uses

func (m *Dense) Dims() (r, c int)

Dims returns the number of rows and columns in the matrix.

func (*Dense) DivElem Uses

func (m *Dense) DivElem(a, b Matrix)

DivElem performs element-wise division of a by b, placing the result in the receiver. DivElem will panic if the two matrices do not have the same shape.

Code:

// Initialize two matrices, a and b.
a := mat.NewDense(2, 2, []float64{
    5, 10,
    15, 20,
})
b := mat.NewDense(2, 2, []float64{
    5, 5,
    5, 5,
})

// Divide the elements of a by b, placing the result into a.
a.DivElem(a, b)

// Print the result using the formatter.
fa := mat.Formatted(a, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("a = %v", fa)

Output:

a = ⎡1  2⎤
    ⎣3  4⎦

func (*Dense) EigenvectorsSym Uses

func (m *Dense) EigenvectorsSym(e *EigenSym)

EigenvectorsSym extracts the eigenvectors of the factorized matrix and stores them in the receiver. Each eigenvector is a column corresponding to the respective eigenvalue returned by e.Values.

EigenvectorsSym panics if the factorization was not successful or if the decomposition did not compute the eigenvectors.

func (*Dense) Exp Uses

func (m *Dense) Exp(a Matrix)

Exp calculates the exponential of the matrix a, e^a, placing the result in the receiver. Exp will panic with matrix.ErrShape if a is not square.

Exp uses the scaling and squaring method described in section 3 of http://www.cs.cornell.edu/cv/researchpdf/19ways+.pdf.

Code:

// Initialize a matrix a with some data.
a := mat.NewDense(2, 2, []float64{
    1, 0,
    0, 1,
})

// Take the exponential of the matrix and place the result in m.
var m mat.Dense
m.Exp(a)

// Print the result using the formatter.
fm := mat.Formatted(&m, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("m = %4.2f", fm)

Output:

m = ⎡2.72  0.00⎤
    ⎣0.00  2.72⎦

func (*Dense) Grow Uses

func (m *Dense) Grow(r, c int) Matrix

Grow returns the receiver expanded by r rows and c columns. If the dimensions of the expanded matrix are outside the capacities of the receiver a new allocation is made, otherwise not. Note the receiver itself is not modified during the call to Grow.

func (*Dense) Inverse Uses

func (m *Dense) Inverse(a Matrix) error

Inverse computes the inverse of the matrix a, storing the result into the receiver. If a is ill-conditioned, a Condition error will be returned. Note that matrix inversion is numerically unstable, and should generally be avoided where possible, for example by using the Solve routines.

Code:

// Initialize two matrices, a and ia.
a := mat.NewDense(2, 2, []float64{
    4, 0,
    0, 4,
})
var ia mat.Dense

// Take the inverse of a and place the result in ia.
ia.Inverse(a)

// Print the result using the formatter.
fa := mat.Formatted(&ia, mat.Prefix("     "), mat.Squeeze())
fmt.Printf("ia = %.2g\n\n", fa)

// Confirm that A * A^-1 = I
var r mat.Dense
r.Mul(a, &ia)
fr := mat.Formatted(&r, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("r = %v\n\n", fr)

// The Inverse operation, however, is numerically unstable,
// and should typically be avoided.
// For example, a common need is to find x = A^-1 * b.
// In this case, the SolveVec method of VecDense
// (if b is a Vector) or Solve method of Dense (if b is a
// matrix) should used instead of computing the Inverse of A.
b := mat.NewDense(2, 2, []float64{
    2, 0,
    0, 2,
})
var x mat.Dense
x.Solve(a, b)

// Print the result using the formatter.
fx := mat.Formatted(&x, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("x = %v", fx)

Output:

ia = ⎡0.25    -0⎤
     ⎣   0  0.25⎦

r = ⎡1  0⎤
    ⎣0  1⎦

x = ⎡0.5    0⎤
    ⎣  0  0.5⎦

func (*Dense) IsZero Uses

func (m *Dense) IsZero() bool

IsZero returns whether the receiver is zero-sized. Zero-sized matrices can be the receiver for size-restricted operations. Dense matrices can be zeroed using Reset.

func (Dense) MarshalBinary Uses

func (m Dense) MarshalBinary() ([]byte, error)

MarshalBinary encodes the receiver into a binary form and returns the result.

Dense is little-endian encoded as follows:

 0 -  7  number of rows    (int64)
 8 - 15  number of columns (int64)
16 - ..  matrix data elements (float64)
         [0,0] [0,1] ... [0,ncols-1]
         [1,0] [1,1] ... [1,ncols-1]
         ...
         [nrows-1,0] ... [nrows-1,ncols-1]

func (Dense) MarshalBinaryTo Uses

func (m Dense) MarshalBinaryTo(w io.Writer) (int, error)

MarshalBinaryTo encodes the receiver into a binary form and writes it into w. MarshalBinaryTo returns the number of bytes written into w and an error, if any.

See MarshalBinary for the on-disk layout.

func (*Dense) Mul Uses

func (m *Dense) Mul(a, b Matrix)

Mul takes the matrix product of a and b, placing the result in the receiver. If the number of columns in a does not equal the number of rows in b, Mul will panic.

Code:

// Initialize two matrices, a and b.
a := mat.NewDense(2, 2, []float64{
    4, 0,
    0, 4,
})
b := mat.NewDense(2, 3, []float64{
    4, 0, 0,
    0, 0, 4,
})

// Take the matrix product of a and b and place the result in c.
var c mat.Dense
c.Mul(a, b)

// Print the result using the formatter.
fc := mat.Formatted(&c, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("c = %v", fc)

Output:

c = ⎡16  0   0⎤
    ⎣ 0  0  16⎦

func (*Dense) MulElem Uses

func (m *Dense) MulElem(a, b Matrix)

MulElem performs element-wise multiplication of a and b, placing the result in the receiver. MulElem will panic if the two matrices do not have the same shape.

Code:

// Initialize two matrices, a and b.
a := mat.NewDense(2, 2, []float64{
    1, 2,
    3, 4,
})
b := mat.NewDense(2, 2, []float64{
    1, 2,
    3, 4,
})

// Multiply the elements of a and b, placing the result into a.
a.MulElem(a, b)

// Print the result using the formatter.
fa := mat.Formatted(a, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("a = %v", fa)

Output:

a = ⎡1   4⎤
    ⎣9  16⎦

func (*Dense) Outer Uses

func (m *Dense) Outer(alpha float64, x, y *VecDense)

Outer calculates the outer product of x and y, and stores the result in the receiver.

m = alpha * x * y'

In order to update an existing matrix, see RankOne.

func (*Dense) Permutation Uses

func (m *Dense) Permutation(r int, swaps []int)

Permutation constructs an r×r permutation matrix with the given row swaps. A permutation matrix has exactly one element equal to one in each row and column and all other elements equal to zero. swaps[i] specifies the row with which i will be swapped, which is equivalent to the non-zero column of row i.

func (*Dense) Pow Uses

func (m *Dense) Pow(a Matrix, n int)

Pow calculates the integral power of the matrix a to n, placing the result in the receiver. Pow will panic if n is negative or if a is not square.

Code:

// Initialize a matrix with some data.
a := mat.NewDense(2, 2, []float64{
    4, 4,
    4, 4,
})

// Take the second power of matrix a and place the result in m.
var m mat.Dense
m.Pow(a, 2)

// Print the result using the formatter.
fm := mat.Formatted(&m, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("m = %v\n\n", fm)

// Take the zeroth power of matrix a and place the result in n.
// We expect an identity matrix of the same size as matrix a.
var n mat.Dense
n.Pow(a, 0)

// Print the result using the formatter.
fn := mat.Formatted(&n, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("n = %v", fn)

Output:

m = ⎡32  32⎤
    ⎣32  32⎦

n = ⎡1  0⎤
    ⎣0  1⎦

func (*Dense) Product Uses

func (m *Dense) Product(factors ...Matrix)

Product calculates the product of the given factors and places the result in the receiver. The order of multiplication operations is optimized to minimize the number of floating point operations on the basis that all matrix multiplications are general.

func (*Dense) RankOne Uses

func (m *Dense) RankOne(a Matrix, alpha float64, x, y *VecDense)

RankOne performs a rank-one update to the matrix a and stores the result in the receiver. If a is zero, see Outer.

m = a + alpha * x * y'

func (*Dense) RawMatrix Uses

func (m *Dense) RawMatrix() blas64.General

RawMatrix returns the underlying blas64.General used by the receiver. Changes to elements in the receiver following the call will be reflected in returned blas64.General.

func (*Dense) RawRowView Uses

func (m *Dense) RawRowView(i int) []float64

RawRowView returns a slice backed by the same array as backing the receiver.

func (*Dense) Reset Uses

func (m *Dense) Reset()

Reset zeros the dimensions of the matrix so that it can be reused as the receiver of a dimensionally restricted operation.

See the Reseter interface for more information.

func (*Dense) RowView Uses

func (m *Dense) RowView(i int) Vector

RowView returns row i of the matrix data represented as a column vector, backed by the matrix data.

See RowViewer for more information.

func (*Dense) Scale Uses

func (m *Dense) Scale(f float64, a Matrix)

Scale multiplies the elements of a by f, placing the result in the receiver.

See the Scaler interface for more information.

Code:

// Initialize a matrix with some data.
a := mat.NewDense(2, 2, []float64{
    4, 4,
    4, 4,
})

// Scale the matrix by a factor of 0.25 and place the result in m.
var m mat.Dense
m.Scale(0.25, a)

// Print the result using the formatter.
fm := mat.Formatted(&m, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("m = %4.3f", fm)

Output:

m = ⎡1.000  1.000⎤
    ⎣1.000  1.000⎦

func (*Dense) Set Uses

func (m *Dense) Set(i, j int, v float64)

Set sets the element at row i, column j to the value v.

func (*Dense) SetCol Uses

func (m *Dense) SetCol(j int, src []float64)

SetCol sets the values in the specified column of the matrix to the values in src. len(src) must equal the number of rows in the receiver.

func (*Dense) SetRawMatrix Uses

func (m *Dense) SetRawMatrix(b blas64.General)

SetRawMatrix sets the underlying blas64.General used by the receiver. Changes to elements in the receiver following the call will be reflected in b.

func (*Dense) SetRow Uses

func (m *Dense) SetRow(i int, src []float64)

SetRow sets the values in the specified rows of the matrix to the values in src. len(src) must equal the number of columns in the receiver.

func (*Dense) Slice Uses

func (m *Dense) Slice(i, k, j, l int) Matrix

Slice returns a new Matrix that shares backing data with the receiver. The returned matrix starts at {i,j} of the receiver and extends k-i rows and l-j columns. The final row in the resulting matrix is k-1 and the final column is l-1. Slice panics with ErrIndexOutOfRange if the slice is outside the capacity of the receiver.

func (*Dense) Solve Uses

func (m *Dense) Solve(a, b Matrix) error

Solve finds a minimum-norm solution to a system of linear equations defined by the matrices a and b. If A is singular or near-singular, a Condition error is returned. Please see the documentation for Condition for more information.

The minimization problem solved depends on the input parameters:

- if m >= n, find X such that ||A*X - B||_2 is minimized,
- if m < n, find the minimum norm solution of A * X = B.

The solution matrix, X, is stored in-place into the receiver.

func (*Dense) Stack Uses

func (m *Dense) Stack(a, b Matrix)

Stack appends the rows of b onto the rows of a, placing the result into the receiver with b placed in the greater indexed rows. Stack will panic if the two input matrices do not have the same number of columns or the constructed stacked matrix is not the same shape as the receiver.

func (*Dense) Sub Uses

func (m *Dense) Sub(a, b Matrix)

Sub subtracts the matrix b from a, placing the result in the receiver. Sub will panic if the two matrices do not have the same shape.

Code:

// Initialize two matrices, a and b.
a := mat.NewDense(2, 2, []float64{
    1, 1,
    1, 1,
})
b := mat.NewDense(2, 2, []float64{
    1, 0,
    0, 1,
})

// Subtract b from a, placing the result into a.
a.Sub(a, b)

// Print the result using the formatter.
fa := mat.Formatted(a, mat.Prefix("    "), mat.Squeeze())
fmt.Printf("a = %v", fa)

Output:

a = ⎡0  1⎤
    ⎣1  0⎦

func (*Dense) T Uses

func (m *Dense) T() Matrix

T performs an implicit transpose by returning the receiver inside a Transpose.

func (*Dense) UnmarshalBinary Uses

func (m *Dense) UnmarshalBinary(data []byte) error

UnmarshalBinary decodes the binary form into the receiver. It panics if the receiver is a non-zero Dense matrix.

See MarshalBinary for the on-disk layout.

Limited checks on the validity of the binary input are performed:

- matrix.ErrShape is returned if the number of rows or columns is negative,
- an error is returned if the resulting Dense matrix is too
big for the current architecture (e.g. a 16GB matrix written by a
64b application and read back from a 32b application.)

UnmarshalBinary does not limit the size of the unmarshaled matrix, and so it should not be used on untrusted data.

func (*Dense) UnmarshalBinaryFrom Uses

func (m *Dense) UnmarshalBinaryFrom(r io.Reader) (int, error)

UnmarshalBinaryFrom decodes the binary form into the receiver and returns the number of bytes read and an error if any. It panics if the receiver is a non-zero Dense matrix.

See MarshalBinary for the on-disk layout.

Limited checks on the validity of the binary input are performed:

- matrix.ErrShape is returned if the number of rows or columns is negative,
- an error is returned if the resulting Dense matrix is too
big for the current architecture (e.g. a 16GB matrix written by a
64b application and read back from a 32b application.)

UnmarshalBinary does not limit the size of the unmarshaled matrix, and so it should not be used on untrusted data.

type Eigen Uses

type Eigen struct {
    // contains filtered or unexported fields
}

Eigen is a type for creating and using the eigenvalue decomposition of a dense matrix.

func (*Eigen) Factorize Uses

func (e *Eigen) Factorize(a Matrix, left, right bool) (ok bool)

Factorize computes the eigenvalues of the square matrix a, and optionally the eigenvectors.

A right eigenvalue/eigenvector combination is defined by

A * x_r = λ * x_r

where x_r is the column vector called an eigenvector, and λ is the corresponding eigenvector.

Similarly, a left eigenvalue/eigenvector combination is defined by

x_l * A = λ * x_l

The eigenvalues, but not the eigenvectors, are the same for both decompositions.

Typically eigenvectors refer to right eigenvectors.

In all cases, Eigen computes the eigenvalues of the matrix. If right and left are true, then the right and left eigenvectors will be computed, respectively. Eigen panics if the input matrix is not square.

Factorize returns whether the decomposition succeeded. If the decomposition failed, methods that require a successful factorization will panic.

func (*Eigen) LeftVectors Uses

func (e *Eigen) LeftVectors() *Dense

LeftVectors returns the left eigenvectors of the decomposition. LeftVectors will panic if the left eigenvectors were not computed during the factorization. or if the factorization was not successful.

See the documentation in lapack64.Geev for the format of the vectors.

BUG: This signature and behavior will change when issue #308 is resolved.

func (*Eigen) Values Uses

func (e *Eigen) Values(dst []complex128) []complex128

Values extracts the eigenvalues of the factorized matrix. If dst is non-nil, the values are stored in-place into dst. In this case dst must have length n, otherwise Values will panic. If dst is nil, then a new slice will be allocated of the proper length and filed with the eigenvalues.

Values panics if the Eigen decomposition was not successful.

func (*Eigen) Vectors Uses

func (e *Eigen) Vectors() *Dense

Vectors returns the right eigenvectors of the decomposition. Vectors will panic if the right eigenvectors were not computed during the factorization, or if the factorization was not successful.

The returned matrix will contain the right eigenvectors of the decomposition in the columns of the n×n matrix in the same order as their eigenvalues. If the j-th eigenvalue is real, then

u_j = VL[:,j],
v_j = VR[:,j],

and if it is not real, then j and j+1 form a complex conjugate pair and the eigenvectors can be recovered as

u_j     = VL[:,j] + i*VL[:,j+1],
u_{j+1} = VL[:,j] - i*VL[:,j+1],
v_j     = VR[:,j] + i*VR[:,j+1],
v_{j+1} = VR[:,j] - i*VR[:,j+1],

where i is the imaginary unit. The computed eigenvectors are normalized to have Euclidean norm equal to 1 and largest component real.

BUG: This signature and behavior will change when issue #308 is resolved.

type EigenSym Uses

type EigenSym struct {
    // contains filtered or unexported fields
}

EigenSym is a type for creating and manipulating the Eigen decomposition of symmetric matrices.

func (*EigenSym) Factorize Uses

func (e *EigenSym) Factorize(a Symmetric, vectors bool) (ok bool)

Factorize computes the eigenvalue decomposition of the symmetric matrix a. The Eigen decomposition is defined as

A = P * D * P^-1

where D is a diagonal matrix containing the eigenvalues of the matrix, and P is a matrix of the eigenvectors of A. If the vectors input argument is false, the eigenvectors are not computed.

Factorize returns whether the decomposition succeeded. If the decomposition failed, methods that require a successful factorization will panic.

func (*EigenSym) Values Uses

func (e *EigenSym) Values(dst []float64) []float64

Values extracts the eigenvalues of the factorized matrix. If dst is non-nil, the values are stored in-place into dst. In this case dst must have length n, otherwise Values will panic. If dst is nil, then a new slice will be allocated of the proper length and filled with the eigenvalues.

Values panics if the Eigen decomposition was not successful.

type Error Uses

type Error struct {
    // contains filtered or unexported fields
}

Error represents matrix handling errors. These errors can be recovered by Maybe wrappers.

func (Error) Error Uses

func (err Error) Error() string

type ErrorStack Uses

type ErrorStack struct {
    Err error

    // StackTrace is the stack trace
    // recovered by Maybe, MaybeFloat
    // or MaybeComplex.
    StackTrace string
}

ErrorStack represents matrix handling errors that have been recovered by Maybe wrappers.

func (ErrorStack) Error Uses

func (err ErrorStack) Error() string

type FormatOption Uses

type FormatOption func(*formatter)

FormatOption is a functional option for matrix formatting.

func DotByte Uses

func DotByte(b byte) FormatOption

DotByte sets the dot character to b. The dot character is used to replace zero elements if the result is printed with the fmt ' ' verb flag. Without a DotByte option, the default dot character is '.'.

func Excerpt Uses

func Excerpt(m int) FormatOption

Excerpt sets the maximum number of rows and columns to print at the margins of the matrix to m. If m is zero or less all elements are printed.

Code:

// Excerpt allows diagnostic display of very large
// matrices and vectors.

// The big matrix is too large to properly print...
big := mat.NewDense(100, 100, nil)
for i := 0; i < 100; i++ {
    big.Set(i, i, 1)
}

// so only print corner excerpts of the matrix.
fmt.Printf("excerpt big identity matrix: %v\n\n",
    mat.Formatted(big, mat.Prefix(" "), mat.Excerpt(3)))

// The long vector is also too large, ...
long := mat.NewVecDense(100, nil)
for i := 0; i < 100; i++ {
    long.SetVec(i, float64(i))
}

// ... so print end excerpts of the vector,
fmt.Printf("excerpt long column vector: %v\n\n",
    mat.Formatted(long, mat.Prefix(" "), mat.Excerpt(3)))
// or its transpose.
fmt.Printf("excerpt long row vector: %v\n",
    mat.Formatted(long.T(), mat.Prefix(" "), mat.Excerpt(3)))

Output:

excerpt big identity matrix: Dims(100, 100)
 ⎡1  0  0  ...  ...  0  0  0⎤
 ⎢0  1  0            0  0  0⎥
 ⎢0  0  1            0  0  0⎥
  .
  .
  .
 ⎢0  0  0            1  0  0⎥
 ⎢0  0  0            0  1  0⎥
 ⎣0  0  0  ...  ...  0  0  1⎦

excerpt long column vector: Dims(100, 1)
 ⎡ 0⎤
 ⎢ 1⎥
 ⎢ 2⎥
  .
  .
  .
 ⎢97⎥
 ⎢98⎥
 ⎣99⎦

excerpt long row vector: Dims(1, 100)
 [ 0   1   2  ...  ...  97  98  99]

func Prefix Uses

func Prefix(p string) FormatOption

Prefix sets the formatted prefix to the string p. Prefix is a string that is prepended to each line of output.

func Squeeze Uses

func Squeeze() FormatOption

Squeeze sets the printing behaviour to minimise column width for each individual column.

type GSVD Uses

type GSVD struct {
    // contains filtered or unexported fields
}

GSVD is a type for creating and using the Generalized Singular Value Decomposition (GSVD) of a matrix.

The factorization is a linear transformation of the data sets from the given variable×sample spaces to reduced and diagonalized "eigenvariable"×"eigensample" spaces.

Code:

// Perform a GSVD factorization on food production/consumption data for the
// three years 1990, 2000 and 2014, for Africa and Latin America/Caribbean.
//
// See Lee et al. doi:10.1371/journal.pone.0030098 and
// Alter at al. doi:10.1073/pnas.0530258100 for more details.
var gsvd mat.GSVD
ok := gsvd.Factorize(FAO.Africa, FAO.LatinAmericaCaribbean, mat.GSVDU|mat.GSVDV|mat.GSVDQ)
if !ok {
    log.Fatal("GSVD factorization failed")
}

u := gsvd.UTo(nil)
v := gsvd.VTo(nil)

s1 := gsvd.ValuesA(nil)
s2 := gsvd.ValuesB(nil)

fmt.Printf("Africa\n\ts1 = %.4f\n\n\tU = %.4f\n\n",
    s1, mat.Formatted(u, mat.Prefix("\t    "), mat.Excerpt(2)))
fmt.Printf("Latin America/Caribbean\n\ts2 = %.4f\n\n\tV = %.4f\n",
    s2, mat.Formatted(v, mat.Prefix("\t    "), mat.Excerpt(2)))

var q mat.Dense
q.Mul(gsvd.ZeroRTo(nil), gsvd.QTo(nil))
fmt.Printf("\nCommon basis vectors\n\n\tQ^T = %.4f\n",
    mat.Formatted(q.T(), mat.Prefix("\t      ")))

// Calculate the antisymmetric angular distances for each eigenvariable.
fmt.Println("\nSignificance:")
for i := 0; i < 3; i++ {
    fmt.Printf("\teigenvar_%d: %+.4f\n", i, math.Atan(s1[i]/s2[i])-math.Pi/4)
}

Output:

Africa
	s1 = [1.0000 0.9344 0.5118]

	U = Dims(21, 21)
	    ⎡-0.0005   0.0142  ...  ...  -0.0060  -0.0055⎤
	    ⎢-0.0010   0.0019             0.0071   0.0075⎥
	     .
	     .
	     .
	    ⎢-0.0007  -0.0024             0.9999  -0.0001⎥
	    ⎣-0.0010  -0.0016  ...  ...  -0.0001   0.9999⎦

Latin America/Caribbean
	s2 = [0.0047 0.3563 0.8591]

	V = Dims(14, 14)
	    ⎡ 0.1362   0.0008  ...  ...   0.0700   0.2636⎤
	    ⎢ 0.1830  -0.0040             0.2908   0.7834⎥
	     .
	     .
	     .
	    ⎢-0.2598  -0.0324             0.9339  -0.2170⎥
	    ⎣-0.8386   0.1494  ...  ...  -0.1639   0.4121⎦

Common basis vectors

	Q^T = ⎡ -8172.4084   -4524.2933    4813.9616⎤
	      ⎢ 22581.8020   12397.1070  -16364.8933⎥
	      ⎣ -8910.8462  -10902.1488   15762.8719⎦

Significance:
	eigenvar_0: +0.7807
	eigenvar_1: +0.4211
	eigenvar_2: -0.2482

func (*GSVD) Factorize Uses

func (gsvd *GSVD) Factorize(a, b Matrix, kind GSVDKind) (ok bool)

Factorize computes the generalized singular value decomposition (GSVD) of the input the r×c matrix A and the p×c matrix B. The singular values of A and B are computed in all cases, while the singular vectors are optionally computed depending on the input kind.

The full singular value decomposition (kind == GSVDU|GSVDV|GSVDQ) deconstructs A and B as

A = U * Σ₁ * [ 0 R ] * Q^T

B = V * Σ₂ * [ 0 R ] * Q^T

where Σ₁ and Σ₂ are r×(k+l) and p×(k+l) diagonal matrices of singular values, and U, V and Q are r×r, p×p and c×c orthogonal matrices of singular vectors. k+l is the effective numerical rank of the matrix [ A^T B^T ]^T.

It is frequently not necessary to compute the full GSVD. Computation time and storage costs can be reduced using the appropriate kind. Either only the singular values can be computed (kind == SVDNone), or in conjunction with specific singular vectors (kind bit set according to matrix.GSVDU, matrix.GSVDV and matrix.GSVDQ).

Factorize returns whether the decomposition succeeded. If the decomposition failed, routines that require a successful factorization will panic.

func (*GSVD) GeneralizedValues Uses

func (gsvd *GSVD) GeneralizedValues(v []float64) []float64

GeneralizedValues returns the generalized singular values of the factorized matrices. If the input slice is non-nil, the values will be stored in-place into the slice. In this case, the slice must have length min(r,c)-k, and GeneralizedValues will panic with matrix.ErrSliceLengthMismatch otherwise. If the input slice is nil, a new slice of the appropriate length will be allocated and returned.

GeneralizedValues will panic if the receiver does not contain a successful factorization.

func (*GSVD) Kind Uses

func (gsvd *GSVD) Kind() GSVDKind

Kind returns the matrix.GSVDKind of the decomposition. If no decomposition has been computed, Kind returns 0.

func (*GSVD) QTo Uses

func (gsvd *GSVD) QTo(dst *Dense) *Dense

QTo extracts the matrix Q from the singular value decomposition, storing the result in-place into dst. Q is size c×c. If dst is nil, a new matrix is allocated. The resulting Q matrix is returned.

QTo will panic if the receiver does not contain a successful factorization.

func (*GSVD) Rank Uses

func (gsvd *GSVD) Rank() (k, l int)

Rank returns the k and l terms of the rank of [ A^T B^T ]^T.

func (*GSVD) SigmaATo Uses

func (gsvd *GSVD) SigmaATo(dst *Dense) *Dense

SigmaATo extracts the matrix Σ₁ from the singular value decomposition, storing the result in-place into dst. Σ₁ is size r×(k+l). If dst is nil, a new matrix is allocated. The resulting SigmaA matrix is returned.

SigmaATo will panic if the receiver does not contain a successful factorization.

func (*GSVD) SigmaBTo Uses

func (gsvd *GSVD) SigmaBTo(dst *Dense) *Dense

SigmaBTo extracts the matrix Σ₂ from the singular value decomposition, storing the result in-place into dst. Σ₂ is size p×(k+l). If dst is nil, a new matrix is allocated. The resulting SigmaB matrix is returned.

SigmaBTo will panic if the receiver does not contain a successful factorization.

func (*GSVD) UTo Uses

func (gsvd *GSVD) UTo(dst *Dense) *Dense

UTo extracts the matrix U from the singular value decomposition, storing the result in-place into dst. U is size r×r. If dst is nil, a new matrix is allocated. The resulting U matrix is returned.

UTo will panic if the receiver does not contain a successful factorization.

func (*GSVD) VTo Uses

func (gsvd *GSVD) VTo(dst *Dense) *Dense

VTo extracts the matrix V from the singular value decomposition, storing the result in-place into dst. V is size p×p. If dst is nil, a new matrix is allocated. The resulting V matrix is returned.

VTo will panic if the receiver does not contain a successful factorization.

func (*GSVD) ValuesA Uses

func (gsvd *GSVD) ValuesA(s []float64) []float64

ValuesA returns the singular values of the factorized A matrix. If the input slice is non-nil, the values will be stored in-place into the slice. In this case, the slice must have length min(r,c)-k, and ValuesA will panic with matrix.ErrSliceLengthMismatch otherwise. If the input slice is nil, a new slice of the appropriate length will be allocated and returned.

ValuesA will panic if the receiver does not contain a successful factorization.

func (*GSVD) ValuesB Uses

func (gsvd *GSVD) ValuesB(s []float64) []float64

ValuesB returns the singular values of the factorized B matrix. If the input slice is non-nil, the values will be stored in-place into the slice. In this case, the slice must have length min(r,c)-k, and ValuesB will panic with matrix.ErrSliceLengthMismatch otherwise. If the input slice is nil, a new slice of the appropriate length will be allocated and returned.

ValuesB will panic if the receiver does not contain a successful factorization.

func (*GSVD) ZeroRTo Uses

func (gsvd *GSVD) ZeroRTo(dst *Dense) *Dense

ZeroRTo extracts the matrix [ 0 R ] from the singular value decomposition, storing the result in-place into dst. [ 0 R ] is size (k+l)×c. If dst is nil, a new matrix is allocated. The resulting ZeroR matrix is returned.

ZeroRTo will panic if the receiver does not contain a successful factorization.

type GSVDKind Uses

type GSVDKind int

GSVDKind specifies the treatment of singular vectors during a GSVD factorization.

const (
    // GSVDU specifies that the U singular vectors should be computed during
    // the decomposition.
    GSVDU GSVDKind = 1 << iota
    // GSVDV specifies that the V singular vectors should be computed during
    // the decomposition.
    GSVDV
    // GSVDQ specifies that the Q singular vectors should be computed during
    // the decomposition.
    GSVDQ

    // GSVDNone specifies that no singular vector should be computed during
    // the decomposition.
    GSVDNone
)

type Grower Uses

type Grower interface {
    Caps() (r, c int)
    Grow(r, c int) Matrix
}

A Grower can grow the size of the represented matrix by the given number of rows and columns. Growing beyond the size given by the Caps method will result in the allocation of a new matrix and copying of the elements. If Grow is called with negative increments it will panic with ErrIndexOutOfRange.

type HOGSVD Uses

type HOGSVD struct {
    // contains filtered or unexported fields
}

HOGSVD is a type for creating and using the Higher Order Generalized Singular Value Decomposition (HOGSVD) of a set of matrices.

The factorization is a linear transformation of the data sets from the given variable×sample spaces to reduced and diagonalized "eigenvariable"×"eigensample" spaces.

Code:

// Perform an HOGSVD factorization on food production/consumption data for the
// three years 1990, 2000 and 2014.
//
// See Ponnapalli et al. doi:10.1371/journal.pone.0028072 and
// Alter at al. doi:10.1073/pnas.0530258100 for more details.
var gsvd mat.HOGSVD
ok := gsvd.Factorize(FAO.Africa, FAO.Asia, FAO.LatinAmericaCaribbean, FAO.Oceania)
if !ok {
    log.Fatalf("HOGSVD factorization failed: %v", gsvd.Err())
}

for i, n := range []string{"Africa", "Asia", "Latin America/Caribbean", "Oceania"} {
    u := gsvd.UTo(nil, i)
    s := gsvd.Values(nil, i)
    fmt.Printf("%s\n\ts_%d = %.4f\n\n\tU_%[2]d = %.4[4]f\n",
        n, i, s, mat.Formatted(u, mat.Prefix("\t      ")))
}

v := gsvd.VTo(nil)
fmt.Printf("\nCommon basis vectors\n\n\tV^T = %.4f",
    mat.Formatted(v.T(), mat.Prefix("\t      ")))

Output:

Africa
	s_0 = [45507.3278 18541.9293 21503.0778]

	U_0 = ⎡-0.0005  -0.0039  -0.0019⎤
	      ⎢-0.0010  -0.0007  -0.0012⎥
	      ⎢-1.0000  -0.0507  -0.9964⎥
	      ⎢-0.0022  -0.2906  -0.0415⎥
	      ⎢ 0.0001  -0.0127  -0.0016⎥
	      ⎢ 0.0003  -0.0067  -0.0010⎥
	      ⎢ 0.0003  -0.0022  -0.0003⎥
	      ⎢-0.0086  -0.9550   0.0734⎥
	      ⎢ 0.0017   0.0002   0.0059⎥
	      ⎢-0.0002  -0.0088  -0.0014⎥
	      ⎢-0.0006  -0.0078  -0.0001⎥
	      ⎢-0.0005  -0.0076   0.0003⎥
	      ⎢ 0.0001  -0.0090   0.0008⎥
	      ⎢-0.0005  -0.0050   0.0029⎥
	      ⎢-0.0011  -0.0078  -0.0012⎥
	      ⎢-0.0014  -0.0058  -0.0002⎥
	      ⎢ 0.0007  -0.0095   0.0020⎥
	      ⎢-0.0008  -0.0081  -0.0009⎥
	      ⎢ 0.0004  -0.0092   0.0006⎥
	      ⎢-0.0007  -0.0079  -0.0006⎥
	      ⎣-0.0011  -0.0076  -0.0010⎦
Asia
	s_1 = [77228.2804 8413.7024 14711.1879]

	U_1 = ⎡ 0.0005  -0.0080   0.0011⎤
	      ⎢ 0.0008  -0.0108   0.0016⎥
	      ⎢-0.9998   0.0612   0.9949⎥
	      ⎢ 0.0007  -0.5734  -0.0468⎥
	      ⎢ 0.0001  -0.0265  -0.0022⎥
	      ⎢ 0.0001  -0.0165  -0.0019⎥
	      ⎢ 0.0000  -0.0070  -0.0013⎥
	      ⎢ 0.0196  -0.8148   0.0893⎥
	      ⎢ 0.0002  -0.0063   0.0012⎥
	      ⎢-0.0001  -0.0135  -0.0013⎥
	      ⎢-0.0004  -0.0135   0.0019⎥
	      ⎢-0.0005  -0.0132   0.0014⎥
	      ⎢ 0.0003  -0.0155   0.0045⎥
	      ⎢-0.0003  -0.0130   0.0025⎥
	      ⎢-0.0007  -0.0105   0.0016⎥
	      ⎢-0.0006  -0.0129   0.0007⎥
	      ⎢-0.0006  -0.0178  -0.0023⎥
	      ⎢-0.0003  -0.0149   0.0016⎥
	      ⎢-0.0001  -0.0134   0.0030⎥
	      ⎢-0.0004  -0.0154   0.0010⎥
	      ⎣-0.0009  -0.0147  -0.0019⎦
Latin America/Caribbean
	s_2 = [274.1364 20736.3116 729.6947]

	U_2 = ⎡ 0.1060  -0.0021   0.0174⎤
	      ⎢ 0.1415  -0.0016   0.0289⎥
	      ⎢ 0.2350  -0.2669  -0.9212⎥
	      ⎢ 0.0290  -0.0118  -0.0429⎥
	      ⎢ 0.0226  -0.0043  -0.0213⎥
	      ⎢ 0.0117  -0.0016  -0.0197⎥
	      ⎢-0.6263  -0.9635   0.2234⎥
	      ⎢ 0.2334  -0.0013   0.1275⎥
	      ⎢-0.0358  -0.0085  -0.0498⎥
	      ⎢-0.1238  -0.0054   0.0313⎥
	      ⎢-0.0421  -0.0059   0.0528⎥
	      ⎢-0.1471  -0.0056   0.0350⎥
	      ⎢-0.2158  -0.0052  -0.0044⎥
	      ⎣-0.6154  -0.0078  -0.2717⎦
Oceania
	s_3 = [8954.1914 6942.6316 17233.0561]

	U_3 = ⎡-0.0080  -0.0012  -0.0040⎤
	      ⎢ 0.0004  -0.0014   0.0001⎥
	      ⎢ 0.9973  -0.0315   0.9991⎥
	      ⎢ 0.0473  -0.7426  -0.0359⎥
	      ⎢ 0.0018  -0.0342  -0.0020⎥
	      ⎢-0.0005  -0.0148  -0.0016⎥
	      ⎢-0.0004  -0.0047  -0.0007⎥
	      ⎢-0.0246  -0.6642  -0.0138⎥
	      ⎢ 0.0003  -0.0287  -0.0023⎥
	      ⎢-0.0011  -0.0148  -0.0014⎥
	      ⎢-0.0108  -0.0198  -0.0039⎥
	      ⎢-0.0149  -0.0183  -0.0048⎥
	      ⎢-0.0178  -0.0208  -0.0075⎥
	      ⎢-0.0266  -0.0063  -0.0016⎥
	      ⎢-0.0012  -0.0234  -0.0006⎥
	      ⎢-0.0084  -0.0184  -0.0030⎥
	      ⎢-0.0232  -0.0191  -0.0124⎥
	      ⎢-0.0072  -0.0226  -0.0035⎥
	      ⎢-0.0150  -0.0144  -0.0045⎥
	      ⎢-0.0068  -0.0227  -0.0034⎥
	      ⎣-0.0127  -0.0136  -0.0049⎦

Common basis vectors

	V^T = ⎡-0.0897  -0.4460  -0.8905⎤
	      ⎢-0.4911  -0.5432  -0.6810⎥
	      ⎣ 0.0644   0.2841   0.9566⎦

func (*HOGSVD) Err Uses

func (gsvd *HOGSVD) Err() error

Err returns the reason for a factorization failure.

func (*HOGSVD) Factorize Uses

func (gsvd *HOGSVD) Factorize(m ...Matrix) (ok bool)

Factorize computes the higher order generalized singular value decomposition (HOGSVD) of the n input r_i×c column tall matrices in m. HOGSV extends the GSVD case from 2 to n input matrices.

M_0 = U_0 * Σ_0 * V^T
M_1 = U_1 * Σ_1 * V^T
.
.
.
M_{n-1} = U_{n-1} * Σ_{n-1} * V^T

where U_i are r_i×c matrices of singular vectors, Σ are c×c matrices singular values, and V is a c×c matrix of singular vectors.

Factorize returns whether the decomposition succeeded. If the decomposition failed, routines that require a successful factorization will panic.

func (*HOGSVD) Len Uses

func (gsvd *HOGSVD) Len() int

Len returns the number of matrices that have been factorized. If Len returns zero, the factorization was not successful.

func (*HOGSVD) UTo Uses

func (gsvd *HOGSVD) UTo(dst *Dense, n int) *Dense

UTo extracts the matrix U_n from the singular value decomposition, storing the result in-place into dst. U_n is size r×c. If dst is nil, a new matrix is allocated. The resulting U matrix is returned.

UTo will panic if the receiver does not contain a successful factorization.

func (*HOGSVD) VTo Uses

func (gsvd *HOGSVD) VTo(dst *Dense) *Dense

VTo extracts the matrix V from the singular value decomposition, storing the result in-place into dst. V is size c×c. If dst is nil, a new matrix is allocated. The resulting V matrix is returned.

VTo will panic if the receiver does not contain a successful factorization.

func (*HOGSVD) Values Uses

func (gsvd *HOGSVD) Values(s []float64, n int) []float64

Values returns the nth set of singular values of the factorized system. If the input slice is non-nil, the values will be stored in-place into the slice. In this case, the slice must have length c, and Values will panic with matrix.ErrSliceLengthMismatch otherwise. If the input slice is nil, a new slice of the appropriate length will be allocated and returned.

Values will panic if the receiver does not contain a successful factorization.

type LQ Uses

type LQ struct {
    // contains filtered or unexported fields
}

LQ is a type for creating and using the LQ factorization of a matrix.

func (*LQ) Cond Uses

func (lq *LQ) Cond() float64

Cond returns the condition number for the factorized matrix. Cond will panic if the receiver does not contain a successful factorization.

func (*LQ) Factorize Uses

func (lq *LQ) Factorize(a Matrix)

Factorize computes the LQ factorization of an m×n matrix a where n <= m. The LQ factorization always exists even if A is singular.

The LQ decomposition is a factorization of the matrix A such that A = L * Q. The matrix Q is an orthonormal n×n matrix, and L is an m×n upper triangular matrix. L and Q can be extracted from the LTo and QTo methods.

func (*LQ) LTo Uses

func (lq *LQ) LTo(dst *Dense) *Dense

LTo extracts the m×n lower trapezoidal matrix from a LQ decomposition. If dst is nil, a new matrix is allocated. The resulting L matrix is returned.

func (*LQ) QTo Uses

func (lq *LQ) QTo(dst *Dense) *Dense

QTo extracts the n×n orthonormal matrix Q from an LQ decomposition. If dst is nil, a new matrix is allocated. The resulting Q matrix is returned.

func (*LQ) Solve Uses

func (lq *LQ) Solve(m *Dense, trans bool, b Matrix) error

Solve finds a minimum-norm solution to a system of linear equations defined by the matrices A and b, where A is an m×n matrix represented in its LQ factorized form. If A is singular or near-singular a Condition error is returned. Please see the documentation for Condition for more information.

The minimization problem solved depends on the input parameters.

If trans == false, find the minimum norm solution of A * X = b.
If trans == true, find X such that ||A*X - b||_2 is minimized.

The solution matrix, X, is stored in place into m.

func (*LQ) SolveVec Uses

func (lq *LQ) SolveVec(v *VecDense, trans bool, b *VecDense) error

SolveVec finds a minimum-norm solution to a system of linear equations. Please see LQ.Solve for the full documentation.

type LU Uses

type LU struct {
    // contains filtered or unexported fields
}

LU is a type for creating and using the LU factorization of a matrix.

func (*LU) Cond Uses

func (lu *LU) Cond() float64

Cond returns the condition number for the factorized matrix. Cond will panic if the receiver does not contain a successful factorization.

func (*LU) Det Uses

func (lu *LU) Det() float64

Det returns the determinant of the matrix that has been factorized. In many expressions, using LogDet will be more numerically stable.

func (*LU) Factorize Uses

func (lu *LU) Factorize(a Matrix)

Factorize computes the LU factorization of the square matrix a and stores the result. The LU decomposition will complete regardless of the singularity of a.

The LU factorization is computed with pivoting, and so really the decomposition is a PLU decomposition where P is a permutation matrix. The individual matrix factors can be extracted from the factorization using the Permutation method on Dense, and the LU LTo and UTo methods.

func (*LU) LTo Uses

func (lu *LU) LTo(dst *TriDense) *TriDense

LTo extracts the lower triangular matrix from an LU factorization. If dst is nil, a new matrix is allocated. The resulting L matrix is returned.

func (*LU) LogDet Uses

func (lu *LU) LogDet() (det float64, sign float64)

LogDet returns the log of the determinant and the sign of the determinant for the matrix that has been factorized. Numerical stability in product and division expressions is generally improved by working in log space.

func (*LU) Pivot Uses

func (lu *LU) Pivot(swaps []int) []int

Pivot returns pivot indices that enable the construction of the permutation matrix P (see Dense.Permutation). If swaps == nil, then new memory will be allocated, otherwise the length of the input must be equal to the size of the factorized matrix.

func (*LU) RankOne Uses

func (lu *LU) RankOne(orig *LU, alpha float64, x, y *VecDense)

RankOne updates an LU factorization as if a rank-one update had been applied to the original matrix A, storing the result into the receiver. That is, if in the original LU decomposition P * L * U = A, in the updated decomposition P * L * U = A + alpha * x * y^T.

func (*LU) Reset Uses

func (lu *LU) Reset()

Reset resets the factorization so that it can be reused as the receiver of a dimensionally restricted operation.

func (*LU) Solve Uses

func (lu *LU) Solve(m *Dense, trans bool, b Matrix) error

Solve solves a system of linear equations using the LU decomposition of a matrix. It computes

A * x = b if trans == false
A^T * x = b if trans == true

In both cases, A is represented in LU factorized form, and the matrix x is stored into m.

If A is singular or near-singular a Condition error is returned. Please see the documentation for Condition for more information.

func (*LU) SolveVec Uses

func (lu *LU) SolveVec(v *VecDense, trans bool, b *VecDense) error

SolveVec solves a system of linear equations using the LU decomposition of a matrix. It computes

A * x = b if trans == false
A^T * x = b if trans == true

In both cases, A is represented in LU factorized form, and the matrix x is stored into v.

If A is singular or near-singular a Condition error is returned. Please see the documentation for Condition for more information.

func (*LU) UTo Uses

func (lu *LU) UTo(dst *TriDense) *TriDense

UTo extracts the upper triangular matrix from an LU factorization. If dst is nil, a new matrix is allocated. The resulting U matrix is returned.

type Matrix Uses

type Matrix interface {
    // Dims returns the dimensions of a Matrix.
    Dims() (r, c int)

    // At returns the value of a matrix element at row i, column j.
    // It will panic if i or j are out of bounds for the matrix.
    At(i, j int) float64

    // T returns the transpose of the Matrix. Whether T returns a copy of the
    // underlying data is implementation dependent.
    // This method may be implemented using the Transpose type, which
    // provides an implicit matrix transpose.
    T() Matrix
}

Matrix is the basic matrix interface type.

type Mutable Uses

type Mutable interface {
    // Set alters the matrix element at row i, column j to v.
    // It will panic if i or j are out of bounds for the matrix.
    Set(i, j int, v float64)

    Matrix
}

Mutable is a matrix interface type that allows elements to be altered.

type MutableBanded Uses

type MutableBanded interface {
    Banded
    SetBand(i, j int, v float64)
}

A MutableBanded can set elements of a band matrix.

type MutableSymBanded Uses

type MutableSymBanded interface {
    Symmetric
    Bandwidth() (kl, ku int)
    SetSymBand(i, j int, v float64)
}

MutableSymBanded is a symmetric band matrix interface type that allows elements to be altered.

type MutableSymmetric Uses

type MutableSymmetric interface {
    Symmetric
    SetSym(i, j int, v float64)
}

A MutableSymmetric can set elements of a symmetric matrix.

type MutableTriangular Uses

type MutableTriangular interface {
    Triangular
    SetTri(i, j int, v float64)
}

A MutableTriangular can set elements of a triangular matrix.

type NonZeroDoer Uses

type NonZeroDoer interface {
    DoNonZero(func(i, j int, v float64))
}

A NonZeroDoer can call a function for each non-zero element of the receiver. The parameters of the function are the element indices and its value.

type QR Uses

type QR struct {
    // contains filtered or unexported fields
}

QR is a type for creating and using the QR factorization of a matrix.

func (*QR) Cond Uses

func (qr *QR) Cond() float64

Cond returns the condition number for the factorized matrix. Cond will panic if the receiver does not contain a successful factorization.

func (*QR) Factorize Uses

func (qr *QR) Factorize(a Matrix)

Factorize computes the QR factorization of an m×n matrix a where m >= n. The QR factorization always exists even if A is singular.

The QR decomposition is a factorization of the matrix A such that A = Q * R. The matrix Q is an orthonormal m×m matrix, and R is an m×n upper triangular matrix. Q and R can be extracted using the QTo and RTo methods.

func (*QR) QTo Uses

func (qr *QR) QTo(dst *Dense) *Dense

QTo extracts the m×m orthonormal matrix Q from a QR decomposition. If dst is nil, a new matrix is allocated. The resulting Q matrix is returned.

func (*QR) RTo Uses

func (qr *QR) RTo(dst *Dense) *Dense

RTo extracts the m×n upper trapezoidal matrix from a QR decomposition. If dst is nil, a new matrix is allocated. The resulting dst matrix is returned.

func (*QR) Solve Uses

func (qr *QR) Solve(m *Dense, trans bool, b Matrix) error

Solve finds a minimum-norm solution to a system of linear equations defined by the matrices A and b, where A is an m×n matrix represented in its QR factorized form. If A is singular or near-singular a Condition error is returned. Please see the documentation for Condition for more information.

The minimization problem solved depends on the input parameters.

If trans == false, find X such that ||A*X - b||_2 is minimized.
If trans == true, find the minimum norm solution of A^T * X = b.

The solution matrix, X, is stored in place into m.

func (*QR) SolveVec Uses

func (qr *QR) SolveVec(v *VecDense, trans bool, b *VecDense) error

SolveVec finds a minimum-norm solution to a system of linear equations. Please see QR.Solve for the full documentation.

type RawBander Uses

type RawBander interface {
    RawBand() blas64.Band
}

A RawBander can return a blas64.Band representation of the receiver. Changes to the blas64.Band.Data slice will be reflected in the original matrix, changes to the Rows, Cols, KL, KU and Stride fields will not.

type RawColViewer Uses

type RawColViewer interface {
    RawColView(j int) []float64
}

A RawColViewer can return a slice of float64 reflecting a column that is backed by the matrix data.

type RawMatrixSetter Uses

type RawMatrixSetter interface {
    SetRawMatrix(a blas64.General)
}

A RawMatrixSetter can set the underlying blas64.General used by the receiver. There is no restriction on the shape of the receiver. Changes to the receiver's elements will be reflected in the blas64.General.Data.

type RawMatrixer Uses

type RawMatrixer interface {
    RawMatrix() blas64.General
}

A RawMatrixer can return a blas64.General representation of the receiver. Changes to the blas64.General.Data slice will be reflected in the original matrix, changes to the Rows, Cols and Stride fields will not.

type RawRowViewer Uses

type RawRowViewer interface {
    RawRowView(i int) []float64
}

A RawRowViewer can return a slice of float64 reflecting a row that is backed by the matrix data.

type RawSymBander Uses

type RawSymBander interface {
    RawSymBand() blas64.SymmetricBand
}

A RawSymBander can return a blas64.SymmetricBand representation of the receiver. Changes to the blas64.SymmetricBand.Data slice will be reflected in the original matrix, changes to the N, K, Stride and Uplo fields will not.

type RawSymmetricer Uses

type RawSymmetricer interface {
    RawSymmetric() blas64.Symmetric
}

A RawSymmetricer can return a view of itself as a BLAS Symmetric matrix.

type RawTriangular Uses

type RawTriangular interface {
    RawTriangular() blas64.Triangular
}

A RawTriangular can return a view of itself as a BLAS Triangular matrix.

type RawVectorer Uses

type RawVectorer interface {
    RawVector() blas64.Vector
}

A RawVectorer can return a blas64.Vector representation of the receiver. Changes to the blas64.Vector.Data slice will be reflected in the original matrix, changes to the Inc field will not.

type Reseter Uses

type Reseter interface {
    Reset()
}

A Reseter can reset the matrix so that it can be reused as the receiver of a dimensionally restricted operation. This is commonly used when the matrix is being used as a workspace or temporary matrix.

If the matrix is a view, using the reset matrix may result in data corruption in elements outside the view.

type RowNonZeroDoer Uses

type RowNonZeroDoer interface {
    DoRowNonZero(i int, fn func(i, j int, v float64))
}

A RowNonZeroDoer can call a function for each non-zero element of a row of the receiver. The parameters of the function are the element indices and its value.

type RowViewer Uses

type RowViewer interface {
    RowView(i int) Vector
}

A RowViewer can return a Vector reflecting a row that is backed by the matrix data. The Vector returned will have length equal to the number of columns.

type SVD Uses

type SVD struct {
    // contains filtered or unexported fields
}

SVD is a type for creating and using the Singular Value Decomposition (SVD) of a matrix.

func (*SVD) Cond Uses

func (svd *SVD) Cond() float64

Cond returns the 2-norm condition number for the factorized matrix. Cond will panic if the receiver does not contain a successful factorization.

func (*SVD) Factorize Uses

func (svd *SVD) Factorize(a Matrix, kind SVDKind) (ok bool)

Factorize computes the singular value decomposition (SVD) of the input matrix A. The singular values of A are computed in all cases, while the singular vectors are optionally computed depending on the input kind.

The full singular value decomposition (kind == SVDFull) deconstructs A as

A = U * Σ * V^T

where Σ is an m×n diagonal matrix of singular vectors, U is an m×m unitary matrix of left singular vectors, and V is an n×n matrix of right singular vectors.

It is frequently not necessary to compute the full SVD. Computation time and storage costs can be reduced using the appropriate kind. Only the singular values can be computed (kind == SVDNone), or a "thin" representation of the singular vectors (kind = SVDThin). The thin representation can save a significant amount of memory if m >> n. See the documentation for the lapack.SVDKind values for more information.

Factorize returns whether the decomposition succeeded. If the decomposition failed, routines that require a successful factorization will panic.

func (*SVD) Kind Uses

func (svd *SVD) Kind() SVDKind

Kind returns the matrix.SVDKind of the decomposition. If no decomposition has been computed, Kind returns 0.

func (*SVD) UTo Uses

func (svd *SVD) UTo(dst *Dense) *Dense

UTo extracts the matrix U from the singular value decomposition, storing the result in-place into dst. U is size m×m if svd.Kind() == SVDFull, of size m×min(m,n) if svd.Kind() == SVDThin, and UTo panics otherwise.

func (*SVD) VTo Uses

func (svd *SVD) VTo(dst *Dense) *Dense

VTo extracts the matrix V from the singular value decomposition, storing the result in-place into dst. V is size n×n if svd.Kind() == SVDFull, of size n×min(m,n) if svd.Kind() == SVDThin, and VTo panics otherwise.

func (*SVD) Values Uses

func (svd *SVD) Values(s []float64) []float64

Values returns the singular values of the factorized matrix in decreasing order. If the input slice is non-nil, the values will be stored in-place into the slice. In this case, the slice must have length min(m,n), and Values will panic with matrix.ErrSliceLengthMismatch otherwise. If the input slice is nil, a new slice of the appropriate length will be allocated and returned.

Values will panic if the receiver does not contain a successful factorization.

type SVDKind Uses

type SVDKind int

SVDKind specifies the treatment of singular vectors during an SVD factorization.

const (
    // SVDNone specifies that no singular vectors should be computed during
    // the decomposition.
    SVDNone SVDKind = iota + 1
    // SVDThin computes the thin singular vectors, that is, it computes
    //  A = U~ * Σ * V~^T
    // where U~ is of size m×min(m,n), Σ is a diagonal matrix of size min(m,n)×min(m,n)
    // and V~ is of size n×min(m,n).
    SVDThin
    // SVDFull computes the full singular value decomposition,
    //  A = U * Σ * V^T
    // where U is of size m×m, Σ is an m×n diagonal matrix, and V is an n×n matrix.
    SVDFull
)

type SymBandDense Uses

type SymBandDense struct {
    // contains filtered or unexported fields
}

SymBandDense represents a symmetric band matrix in dense storage format.

func NewDiagonal Uses

func NewDiagonal(n int, data []float64) *SymBandDense

NewDiagonal is a convenience function that returns a diagonal matrix represented by a SymBandDense. The length of data must be n or data must be nil, otherwise NewDiagonal will panic.

func NewSymBandDense Uses

func NewSymBandDense(n, k int, data []float64) *SymBandDense

NewSymBandDense creates a new SymBand matrix with n rows and columns. If data == nil, a new slice is allocated for the backing slice. If len(data) == n*(k+1), data is used as the backing slice, and changes to the elements of the returned SymBandDense will be reflected in data. If neither of these is true, NewSymBandDense will panic. k must be at least zero and less than n, otherwise NewBandDense will panic.

The data must be arranged in row-major order constructed by removing the zeros from the rows outside the band and aligning the diagonals. SymBandDense matrices are stored in the upper triangle. For example, the matrix

1  2  3  0  0  0
2  4  5  6  0  0
3  5  7  8  9  0
0  6  8 10 11 12
0  0  9 11 13 14
0  0  0 12 14 15

becomes (* entries are never accessed)

 1  2  3
 4  5  6
 7  8  9
10 11 12
13 14  *
15  *  *

which is passed to NewBandDense as []float64{1, 2, 3, 4, ...} with k=2. Only the values in the band portion of the matrix are used.

func (*SymBandDense) At Uses

func (s *SymBandDense) At(i, j int) float64

At returns the element at row i, column j.

func (*SymBandDense) Bandwidth Uses

func (s *SymBandDense) Bandwidth() (kl, ku int)

Bandwidth returns the bandwidths of the matrix.

func (*SymBandDense) Dims Uses

func (s *SymBandDense) Dims() (r, c int)

Dims returns the number of rows and columns in the matrix.

func (*SymBandDense) DoColNonZero Uses

func (s *SymBandDense) DoColNonZero(j int, fn func(i, j int, v float64))

DoColNonZero calls the function fn for each of the non-zero elements of column j of s. The function fn takes a row/column index and the element value of s at (i, j).

func (*SymBandDense) DoNonZero Uses

func (s *SymBandDense) DoNonZero(fn func(i, j int, v float64))

DoNonZero calls the function fn for each of the non-zero elements of s. The function fn takes a row/column index and the element value of s at (i, j).

func (*SymBandDense) DoRowNonZero Uses

func (s *SymBandDense) DoRowNonZero(i int, fn func(i, j int, v float64))

DoRowNonZero calls the function fn for each of the non-zero elements of row i of s. The function fn takes a row/column index and the element value of s at (i, j).

func (*SymBandDense) RawSymBand Uses

func (s *SymBandDense) RawSymBand() blas64.SymmetricBand

RawSymBand returns the underlying blas64.SymBand used by the receiver. Changes to elements in the receiver following the call will be reflected in returned blas64.SymBand.

func (*SymBandDense) SetSymBand Uses

func (s *SymBandDense) SetSymBand(i, j int, v float64)

SetSymBand sets the element at row i, column j to the value v. It panics if the location is outside the appropriate region of the matrix.

func (*SymBandDense) Symmetric Uses

func (s *SymBandDense) Symmetric() int

Symmetric returns the size of the receiver.

func (*SymBandDense) T Uses

func (s *SymBandDense) T() Matrix

T implements the Matrix interface. Symmetric matrices, by definition, are equal to their transpose, and this is a no-op.

func (*SymBandDense) TBand Uses

func (s *SymBandDense) TBand() Banded

TBand implements the Banded interface.

type SymDense Uses

type SymDense struct {
    // contains filtered or unexported fields
}

SymDense is a symmetric matrix that uses dense storage. SymDense matrices are stored in the upper triangle.

func NewSymDense Uses

func NewSymDense(n int, data []float64) *SymDense

NewSymDense creates a new Symmetric matrix with n rows and columns. If data == nil, a new slice is allocated for the backing slice. If len(data) == n*n, data is used as the backing slice, and changes to the elements of the returned SymDense will be reflected in data. If neither of these is true, NewSymDense will panic.

The data must be arranged in row-major order, i.e. the (i*c + j)-th element in the data slice is the {i, j}-th element in the matrix. Only the values in the upper triangular portion of the matrix are used.

func (*SymDense) AddSym Uses

func (s *SymDense) AddSym(a, b Symmetric)

func (*SymDense) At Uses

func (s *SymDense) At(i, j int) float64

At returns the element at row i and column j.

func (*SymDense) Caps Uses

func (s *SymDense) Caps() (r, c int)

Caps returns the number of rows and columns in the backing matrix.

func (*SymDense) CopySym Uses

func (s *SymDense) CopySym(a Symmetric) int

func (*SymDense) Dims Uses

func (s *SymDense) Dims() (r, c int)

Dims returns the number of rows and columns in the matrix.

func (*SymDense) GrowSquare Uses

func (s *SymDense) GrowSquare(n int) Matrix

GrowSquare returns the receiver expanded by n rows and n columns. If the dimensions of the expanded matrix are outside the capacity of the receiver a new allocation is made, otherwise not. Note that the receiver itself is not modified during the call to GrowSquare.

func (*SymDense) IsZero Uses

func (s *SymDense) IsZero() bool

IsZero returns whether the receiver is zero-sized. Zero-sized matrices can be the receiver for size-restricted operations. SymDense matrices can be zeroed using Reset.

func (*SymDense) PowPSD Uses

func (s *SymDense) PowPSD(a Symmetric, pow float64) error

PowPSD computes a^pow where a is a positive symmetric definite matrix.

PowPSD returns an error if the matrix is not not positive symmetric definite or the Eigendecomposition is not successful.

func (*SymDense) RankTwo Uses

func (s *SymDense) RankTwo(a Symmetric, alpha float64, x, y *VecDense)

RankTwo performs a symmmetric rank-two update to the matrix a and stores the result in the receiver

m = a + alpha * (x * y' + y * x')

func (*SymDense) RawSymmetric Uses

func (s *SymDense) RawSymmetric() blas64.Symmetric

RawSymmetric returns the matrix as a blas64.Symmetric. The returned value must be stored in upper triangular format.

func (*SymDense) Reset Uses

func (s *SymDense) Reset()

Reset zeros the dimensions of the matrix so that it can be reused as the receiver of a dimensionally restricted operation.

See the Reseter interface for more information.

func (*SymDense) ScaleSym Uses

func (s *SymDense) ScaleSym(f float64, a Symmetric)

ScaleSym multiplies the elements of a by f, placing the result in the receiver.

func (*SymDense) SetRawSymmetric Uses

func (s *SymDense) SetRawSymmetric(b blas64.Symmetric)

SetRawSymmetric sets the underlying blas64.Symmetric used by the receiver. Changes to elements in the receiver following the call will be reflected in b. SetRawSymmetric will panic if b is not an upper-encoded symmetric matrix.

func (*SymDense) SetSym Uses

func (s *SymDense) SetSym(i, j int, v float64)

SetSym sets the elements at (i,j) and (j,i) to the value v.

func (*SymDense) SliceSquare Uses

func (s *SymDense) SliceSquare(i, k int) Matrix

SliceSquare returns a new Matrix that shares backing data with the receiver. The returned matrix starts at {i,i} of the receiver and extends k-i rows and columns. The final row and column in the resulting matrix is k-1. SliceSquare panics with ErrIndexOutOfRange if the slice is outside the capacity of the receiver.

func (*SymDense) SubsetSym Uses

func (s *SymDense) SubsetSym(a Symmetric, set []int)

SubsetSym extracts a subset of the rows and columns of the matrix a and stores the result in-place into the receiver. The resulting matrix size is len(set)×len(set). Specifically, at the conclusion of SubsetSym, s.At(i, j) equals a.At(set[i], set[j]). Note that the supplied set does not have to be a strict subset, dimension repeats are allowed.

Code:

n := 5
s := mat.NewSymDense(5, nil)
count := 1.0
for i := 0; i < n; i++ {
    for j := i; j < n; j++ {
        s.SetSym(i, j, count)
        count++
    }
}
fmt.Println("Original matrix:")
fmt.Printf("%0.4v\n\n", mat.Formatted(s))

// Take the subset {0, 2, 4}
var sub mat.SymDense
sub.SubsetSym(s, []int{0, 2, 4})
fmt.Println("Subset {0, 2, 4}")
fmt.Printf("%0.4v\n\n", mat.Formatted(&sub))

// Take the subset {0, 0, 4}
sub.SubsetSym(s, []int{0, 0, 4})
fmt.Println("Subset {0, 0, 4}")
fmt.Printf("%0.4v\n\n", mat.Formatted(&sub))

Output:

Original matrix:
⎡ 1   2   3   4   5⎤
⎢ 2   6   7   8   9⎥
⎢ 3   7  10  11  12⎥
⎢ 4   8  11  13  14⎥
⎣ 5   9  12  14  15⎦

Subset {0, 2, 4}
⎡ 1   3   5⎤
⎢ 3  10  12⎥
⎣ 5  12  15⎦

Subset {0, 0, 4}
⎡ 1   1   5⎤
⎢ 1   1   5⎥
⎣ 5   5  15⎦

func (*SymDense) SymOuterK Uses

func (s *SymDense) SymOuterK(alpha float64, x Matrix)

SymOuterK calculates the outer product of x with itself and stores the result into the receiver. It is equivalent to the matrix multiplication

s = alpha * x * x'.

In order to update an existing matrix, see SymRankOne.

func (*SymDense) SymRankK Uses

func (s *SymDense) SymRankK(a Symmetric, alpha float64, x Matrix)

SymRankK performs a symmetric rank-k update to the matrix a and stores the result into the receiver. If a is zero, see SymOuterK.

s = a + alpha * x * x'

func (*SymDense) SymRankOne Uses

func (s *SymDense) SymRankOne(a Symmetric, alpha float64, x *VecDense)

SymRankOne performs a symetric rank-one update to the matrix a and stores the result in the receiver

s = a + alpha * x * x'

func (*SymDense) Symmetric Uses

func (s *SymDense) Symmetric() int

func (*SymDense) T Uses

func (s *SymDense) T() Matrix

T implements the Matrix interface. Symmetric matrices, by definition, are equal to their transpose, and this is a no-op.

type Symmetric Uses

type Symmetric interface {
    Matrix
    // Symmetric returns the number of rows/columns in the matrix.
    Symmetric() int
}

Symmetric represents a symmetric matrix (where the element at {i, j} equals the element at {j, i}). Symmetric matrices are always square.

type Transpose Uses

type Transpose struct {
    Matrix Matrix
}

Transpose is a type for performing an implicit matrix transpose. It implements the Matrix interface, returning values from the transpose of the matrix within.

func (Transpose) At Uses

func (t Transpose) At(i, j int) float64

At returns the value of the element at row i and column j of the transposed matrix, that is, row j and column i of the Matrix field.

func (Transpose) Dims Uses

func (t Transpose) Dims() (r, c int)

Dims returns the dimensions of the transposed matrix. The number of rows returned is the number of columns in the Matrix field, and the number of columns is the number of rows in the Matrix field.

func (Transpose) T Uses

func (t Transpose) T() Matrix

T performs an implicit transpose by returning the Matrix field.

func (Transpose) Untranspose Uses

func (t Transpose) Untranspose() Matrix

Untranspose returns the Matrix field.

type TransposeBand Uses

type TransposeBand struct {
    Banded Banded
}

TransposeBand is a type for performing an implicit transpose of a band matrix. It implements the Banded interface, returning values from the transpose of the matrix within.

func (TransposeBand) At Uses

func (t TransposeBand) At(i, j int) float64

At returns the value of the element at row i and column j of the transposed matrix, that is, row j and column i of the Banded field.

func (TransposeBand) Bandwidth Uses

func (t TransposeBand) Bandwidth() (kl, ku int)

Bandwidth returns the lower and upper bandwidth values for the transposed matrix.

func (TransposeBand) Dims Uses

func (t TransposeBand) Dims() (r, c int)

Dims returns the dimensions of the transposed matrix.

func (TransposeBand) T Uses

func (t TransposeBand) T() Matrix

T performs an implicit transpose by returning the Banded field.

func (TransposeBand) TBand Uses

func (t TransposeBand) TBand() Banded

TBand performs an implicit transpose by returning the Banded field.

func (TransposeBand) Untranspose Uses

func (t TransposeBand) Untranspose() Matrix

Untranspose returns the Banded field.

func (TransposeBand) UntransposeBand Uses

func (t TransposeBand) UntransposeBand() Banded

UntransposeBand returns the Banded field.

type TransposeTri Uses

type TransposeTri struct {
    Triangular Triangular
}

TransposeTri is a type for performing an implicit transpose of a Triangular matrix. It implements the Triangular interface, returning values from the transpose of the matrix within.

func (TransposeTri) At Uses

func (t TransposeTri) At(i, j int) float64

At returns the value of the element at row i and column j of the transposed matrix, that is, row j and column i of the Triangular field.

func (TransposeTri) Dims Uses

func (t TransposeTri) Dims() (r, c int)

Dims returns the dimensions of the transposed matrix. Triangular matrices are square and thus this is the same size as the original Triangular.

func (TransposeTri) T Uses

func (t TransposeTri) T() Matrix

T performs an implicit transpose by returning the Triangular field.

func (TransposeTri) TTri Uses

func (t TransposeTri) TTri() Triangular

TTri performs an implicit transpose by returning the Triangular field.

func (TransposeTri) Triangle Uses

func (t TransposeTri) Triangle() (int, TriKind)

Triangle returns the number of rows/columns in the matrix and its orientation.

func (TransposeTri) Untranspose Uses

func (t TransposeTri) Untranspose() Matrix

Untranspose returns the Triangular field.

func (TransposeTri) UntransposeTri Uses

func (t TransposeTri) UntransposeTri() Triangular

type TriDense Uses

type TriDense struct {
    // contains filtered or unexported fields
}

TriDense represents an upper or lower triangular matrix in dense storage format.

func NewTriDense Uses

func NewTriDense(n int, kind TriKind, data []float64) *TriDense

NewTriDense creates a new Triangular matrix with n rows and columns. If data == nil, a new slice is allocated for the backing slice. If len(data) == n*n, data is used as the backing slice, and changes to the elements of the returned TriDense will be reflected in data. If neither of these is true, NewTriDense will panic.

The data must be arranged in row-major order, i.e. the (i*c + j)-th element in the data slice is the {i, j}-th element in the matrix. Only the values in the triangular portion corresponding to kind are used.

func (*TriDense) At Uses

func (t *TriDense) At(i, j int) float64

At returns the element at row i, column j.

func (*TriDense) Copy Uses

func (t *TriDense) Copy(a Matrix) (r, c int)

Copy makes a copy of elements of a into the receiver. It is similar to the built-in copy; it copies as much as the overlap between the two matrices and returns the number of rows and columns it copied. Only elements within the receiver's non-zero triangle are set.

See the Copier interface for more information.

func (*TriDense) Dims Uses

func (t *TriDense) Dims() (r, c int)

func (*TriDense) DoColNonZero Uses

func (t *TriDense) DoColNonZero(j int, fn func(i, j int, v float64))

DoColNonZero calls the function fn for each of the non-zero elements of column j of t. The function fn takes a row/column index and the element value of t at (i, j).

func (*TriDense) DoNonZero Uses

func (t *TriDense) DoNonZero(fn func(i, j int, v float64))

DoNonZero calls the function fn for each of the non-zero elements of t. The function fn takes a row/column index and the element value of t at (i, j).

func (*TriDense) DoRowNonZero Uses

func (t *TriDense) DoRowNonZero(i int, fn func(i, j int, v float64))

DoRowNonZero calls the function fn for each of the non-zero elements of row i of t. The function fn takes a row/column index and the element value of t at (i, j).

func (*TriDense) InverseTri Uses

func (t *TriDense) InverseTri(a Triangular) error

InverseTri computes the inverse of the triangular matrix a, storing the result into the receiver. If a is ill-conditioned, a Condition error will be returned. Note that matrix inversion is numerically unstable, and should generally be avoided where possible, for example by using the Solve routines.

func (*TriDense) IsZero Uses

func (t *TriDense) IsZero() bool

IsZero returns whether the receiver is zero-sized. Zero-sized matrices can be the receiver for size-restricted operations. TriDense matrices can be zeroed using Reset.

func (*TriDense) MulTri Uses

func (t *TriDense) MulTri(a, b Triangular)

MulTri takes the product of triangular matrices a and b and places the result in the receiver. The size of a and b must match, and they both must have the same TriKind, or Mul will panic.

func (*TriDense) RawTriangular Uses

func (t *TriDense) RawTriangular() blas64.Triangular

func (*TriDense) Reset Uses

func (t *TriDense) Reset()

Reset zeros the dimensions of the matrix so that it can be reused as the receiver of a dimensionally restricted operation.

See the Reseter interface for more information.

func (*TriDense) ScaleTri Uses

func (t *TriDense) ScaleTri(f float64, a Triangular)

ScaleTri multiplies the elements of a by f, placing the result in the receiver. If the receiver is non-zero, the size and kind of the receiver must match the input, or ScaleTri will panic.

func (*TriDense) SetTri Uses

func (t *TriDense) SetTri(i, j int, v float64)

SetTri sets the element at row i, column j to the value v. It panics if the location is outside the appropriate half of the matrix.

func (*TriDense) T Uses

func (t *TriDense) T() Matrix

T performs an implicit transpose by returning the receiver inside a Transpose.

func (*TriDense) TTri Uses

func (t *TriDense) TTri() Triangular

TTri performs an implicit transpose by returning the receiver inside a TransposeTri.

func (*TriDense) Triangle Uses

func (t *TriDense) Triangle() (n int, kind TriKind)

Triangle returns the dimension of t and its orientation. The returned orientation is only valid when n is not zero.

type TriKind Uses

type TriKind bool

TriKind represents the triangularity of the matrix.

const (
    // Upper specifies an upper triangular matrix.
    Upper TriKind = true
    // Lower specifies a lower triangular matrix.
    Lower TriKind = false
)

type Triangular Uses

type Triangular interface {
    Matrix
    // Triangular returns the number of rows/columns in the matrix and its
    // orientation.
    Triangle() (n int, kind TriKind)

    // TTri is the equivalent of the T() method in the Matrix interface but
    // guarantees the transpose is of triangular type.
    TTri() Triangular
}

Triangular represents a triangular matrix. Triangular matrices are always square.

type Unconjugator Uses

type Unconjugator interface {

    // Unconjugate returns the underlying Matrix stored for the implicit
    // conjugate transpose.
    Unconjugate() CMatrix
}

Unconjugator is a type that can undo an implicit conjugate transpose.

type UntransposeBander Uses

type UntransposeBander interface {
    // Untranspose returns the underlying Banded stored for the implicit transpose.
    UntransposeBand() Banded
}

UntransposeBander is a type that can undo an implicit band transpose.

type UntransposeTrier Uses

type UntransposeTrier interface {
    // Untranspose returns the underlying Triangular stored for the implicit transpose.
    UntransposeTri() Triangular
}

UntransposeTrier is a type that can undo an implicit triangular transpose.

type Untransposer Uses

type Untransposer interface {

    // Untranspose returns the underlying Matrix stored for the implicit transpose.
    Untranspose() Matrix
}

Untransposer is a type that can undo an implicit transpose.

type VecDense Uses

type VecDense struct {
    // contains filtered or unexported fields
}

VecDense represents a column vector.

func NewVecDense Uses

func NewVecDense(n int, data []float64) *VecDense

NewVecDense creates a new VecDense of length n. If data == nil, a new slice is allocated for the backing slice. If len(data) == n, data is used as the backing slice, and changes to the elements of the returned VecDense will be reflected in data. If neither of these is true, NewVecDense will panic.

func (*VecDense) AddScaledVec Uses

func (v *VecDense) AddScaledVec(a *VecDense, alpha float64, b *VecDense)

AddScaledVec adds the vectors a and alpha*b, placing the result in the receiver.

func (*VecDense) AddVec Uses

func (v *VecDense) AddVec(a, b *VecDense)

AddVec adds the vectors a and b, placing the result in the receiver.

func (*VecDense) At Uses

func (v *VecDense) At(i, j int) float64

At returns the element at row i. It panics if i is out of bounds or if j is not zero.

func (*VecDense) Cap Uses

func (v *VecDense) Cap() int

Cap returns the capacity of the vector.

func (*VecDense) Caps Uses

func (v *VecDense) Caps() (r, c int)

Caps returns the number of rows and columns in the backing matrix. Columns is always 1 for a non-Reset vector.

func (*VecDense) CloneVec Uses

func (v *VecDense) CloneVec(a *VecDense)

CloneVec makes a copy of a into the receiver, overwriting the previous value of the receiver.

func (*VecDense) ColViewOf Uses

func (v *VecDense) ColViewOf(m RawMatrixer, j int)

ColViewOf reflects the column j of the RawMatrixer m, into the receiver backed by the same underlying data. The length of the receiver must either be zero or match the number of rows in m.

func (*VecDense) CopyVec Uses

func (v *VecDense) CopyVec(a *VecDense) int

CopyVec makes a copy of elements of a into the receiver. It is similar to the built-in copy; it copies as much as the overlap between the two vectors and returns the number of elements it copied.

func (*VecDense) Dims Uses

func (v *VecDense) Dims() (r, c int)

Dims returns the number of rows and columns in the matrix. Columns is always 1 for a non-Reset vector.

func (*VecDense) DivElemVec Uses

func (v *VecDense) DivElemVec(a, b *VecDense)

DivElemVec performs element-wise division of a by b, placing the result in the receiver.

func (*VecDense) IsZero Uses

func (v *VecDense) IsZero() bool

IsZero returns whether the receiver is zero-sized. Zero-sized vectors can be the receiver for size-restricted operations. VecDenses can be zeroed using Reset.

func (*VecDense) Len Uses

func (v *VecDense) Len() int

Len returns the length of the vector.

func (VecDense) MarshalBinary Uses

func (v VecDense) MarshalBinary() ([]byte, error)

MarshalBinary encodes the receiver into a binary form and returns the result.

VecDense is little-endian encoded as follows:

0 -  7  number of elements     (int64)
8 - ..  vector's data elements (float64)

func (VecDense) MarshalBinaryTo Uses

func (v VecDense) MarshalBinaryTo(w io.Writer) (int, error)

MarshalBinaryTo encodes the receiver into a binary form, writes it to w and returns the number of bytes written and an error if any.

See MarshalBainry for the on-disk format.

func (*VecDense) MulElemVec Uses

func (v *VecDense) MulElemVec(a, b *VecDense)

MulElemVec performs element-wise multiplication of a and b, placing the result in the receiver.

func (*VecDense) MulVec Uses

func (v *VecDense) MulVec(a Matrix, b *VecDense)

MulVec computes a * b. The result is stored into the receiver. MulVec panics if the number of columns in a does not equal the number of rows in b.

func (*VecDense) RawVector Uses

func (v *VecDense) RawVector() blas64.Vector

func (*VecDense) Reset Uses

func (v *VecDense) Reset()

Reset zeros the length of the vector so that it can be reused as the receiver of a dimensionally restricted operation.

See the Reseter interface for more information.

func (*VecDense) RowViewOf Uses

func (v *VecDense) RowViewOf(m RawMatrixer, i int)

RowViewOf reflects the row i of the RawMatrixer m, into the receiver backed by the same underlying data. The length of the receiver must either be zero or match the number of columns in m.

func (*VecDense) ScaleVec Uses

func (v *VecDense) ScaleVec(alpha float64, a *VecDense)

ScaleVec scales the vector a by alpha, placing the result in the receiver.

func (*VecDense) SetVec Uses

func (v *VecDense) SetVec(i int, val float64)

SetVec sets the element at row i to the value val. It panics if i is out of bounds.

func (*VecDense) SliceVec Uses

func (v *VecDense) SliceVec(i, k int) *VecDense

SliceVec returns a new VecDense that shares backing data with the receiver. The returned matrix starts at i of the receiver and extends k-i elements. SliceVec panics with ErrIndexOutOfRange if the slice is outside the capacity of the receiver.

func (*VecDense) SolveVec Uses

func (v *VecDense) SolveVec(a Matrix, b *VecDense) error

SolveVec finds a minimum-norm solution to a system of linear equations defined by the matrix a and the right-hand side vector b. If A is singular or near-singular, a Condition error is returned. Please see the documentation for Dense.Solve for more information.

func (*VecDense) SubVec Uses

func (v *VecDense) SubVec(a, b *VecDense)

SubVec subtracts the vector b from a, placing the result in the receiver.

func (*VecDense) T Uses

func (v *VecDense) T() Matrix

T performs an implicit transpose by returning the receiver inside a Transpose.

func (*VecDense) UnmarshalBinary Uses

func (v *VecDense) UnmarshalBinary(data []byte) error

UnmarshalBinary decodes the binary form into the receiver. It panics if the receiver is a non-zero VecDense.

See MarshalBinary for the on-disk layout.

Limited checks on the validity of the binary input are performed:

- matrix.ErrShape is returned if the number of rows is negative,
- an error is returned if the resulting VecDense is too
big for the current architecture (e.g. a 16GB vector written by a
64b application and read back from a 32b application.)

UnmarshalBinary does not limit the size of the unmarshaled vector, and so it should not be used on untrusted data.

func (*VecDense) UnmarshalBinaryFrom Uses

func (v *VecDense) UnmarshalBinaryFrom(r io.Reader) (int, error)

UnmarshalBinaryFrom decodes the binary form into the receiver, from the io.Reader and returns the number of bytes read and an error if any. It panics if the receiver is a non-zero VecDense.

See MarshalBinary for the on-disk layout. See UnmarshalBinary for the list of sanity checks performed on the input.

type Vector Uses

type Vector interface {
    Matrix
    Len() int
}

Vector is a column vector.

Bugs

The computation of the 1-norm and ∞-norm for non-square matrices is innacurate, although is typically the right order of magnitude. See https://github.com/xianyi/OpenBLAS/issues/636. While the value returned will change with the resolution of this bug, the result from Cond will match the condition number used internally.

Package mat imports 15 packages (graph) and is imported by 37 packages. Updated 2017-12-11. Refresh now. Tools for package owners.