Documentation ¶
Index ¶
- Constants
- func GSS(min, max float64, iters int, f func(float64) float64) float64
- type BiCGSTAB
- type BiCGSTABSolver
- type FiniteVector
- type GridSearch2D
- type GridSearch3D
- type KMeans
- type LargeLinearSolver
- type LineSearch
- type Matrix2
- func (m *Matrix2) Add(m1 *Matrix2) *Matrix2
- func (m *Matrix2) Det() float64
- func (m *Matrix2) Eigenvalues() [2]complex128
- func (m *Matrix2) Inverse() *Matrix2
- func (m *Matrix2) InvertInPlace()
- func (m *Matrix2) InvertInPlaceDet(det float64)
- func (m *Matrix2) Mul(m1 *Matrix2) *Matrix2
- func (m *Matrix2) MulColumn(c Vec2) Vec2
- func (m *Matrix2) MulColumnInv(c Vec2, det float64) Vec2
- func (m *Matrix2) SVD(u, s, v *Matrix2)
- func (m *Matrix2) Scale(s float64)
- func (m *Matrix2) Transpose() *Matrix2
- type Matrix3
- func (m *Matrix3) Add(m1 *Matrix3) *Matrix3
- func (m *Matrix3) Cols() [3]Vec3
- func (m *Matrix3) Det() float64
- func (m *Matrix3) Eigenvalues() [3]complex128
- func (m *Matrix3) Inverse() *Matrix3
- func (m *Matrix3) InvertInPlace()
- func (m *Matrix3) InvertInPlaceDet(det float64)
- func (m *Matrix3) Mul(m1 *Matrix3) *Matrix3
- func (m *Matrix3) MulColumn(c Vec3) Vec3
- func (m *Matrix3) MulColumnInv(c Vec3, det float64) Vec3
- func (m *Matrix3) Rows() [3]Vec3
- func (m *Matrix3) SVD(u, s, v *Matrix3)
- func (m *Matrix3) Scale(s float64)
- func (m *Matrix3) Transpose() *Matrix3
- type Matrix4
- func (m *Matrix4) Add(m1 *Matrix4) *Matrix4
- func (m *Matrix4) CharPoly() Polynomial
- func (m *Matrix4) Cols() [4]Vec4
- func (m *Matrix4) Det() float64
- func (m *Matrix4) Mul(m1 *Matrix4) *Matrix4
- func (m *Matrix4) MulColumn(v Vec4) Vec4
- func (m *Matrix4) Rows() [4]Vec4
- func (m *Matrix4) SVD(u, s, v *Matrix4)
- func (m *Matrix4) Scale(s float64) *Matrix4
- func (m *Matrix4) Sub(m1 *Matrix4) *Matrix4
- func (m *Matrix4) Transpose() *Matrix4
- type Polynomial
- func (p Polynomial) Add(p1 Polynomial) Polynomial
- func (p Polynomial) Derivative() Polynomial
- func (p Polynomial) Eval(x float64) float64
- func (p Polynomial) IterRealRoots(f func(x float64) bool)
- func (p Polynomial) Mul(p1 Polynomial) Polynomial
- func (p Polynomial) RealRoots() []float64
- func (p Polynomial) Scale(c float64) Polynomial
- func (p Polynomial) String() string
- type RecursiveLineSearch
- type SparseCholesky
- type SparseMatrix
- func (s *SparseMatrix) Apply(x Vec) Vec
- func (s *SparseMatrix) ApplyVec2(x []Vec2) []Vec2
- func (s *SparseMatrix) ApplyVec3(x []Vec3) []Vec3
- func (s *SparseMatrix) Iterate(row int, f func(col int, x float64))
- func (s *SparseMatrix) Permute(perm []int) *SparseMatrix
- func (s *SparseMatrix) RCM() []int
- func (s *SparseMatrix) Set(row, col int, x float64)
- func (s *SparseMatrix) Transpose() *SparseMatrix
- type Vec
- func (v Vec) Add(v1 Vec) Vec
- func (v Vec) At(i int) float64
- func (v Vec) Dist(v1 Vec) float64
- func (v Vec) DistSquared(v1 Vec) float64
- func (v Vec) Dot(v1 Vec) float64
- func (v Vec) Len() int
- func (v Vec) Norm() float64
- func (v Vec) NormSquared() float64
- func (v Vec) Normalize() Vec
- func (v Vec) ProjectOut(v1 Vec) Vec
- func (v Vec) Scale(s float64) Vec
- func (v Vec) Sub(v1 Vec) Vec
- func (v Vec) WithDim(i int, x float64) Vec
- func (v Vec) Zeros() Vec
- type Vec2
- func (v Vec2) Add(v1 Vec2) Vec2
- func (v Vec2) At(i int) float64
- func (v Vec2) Dist(v1 Vec2) float64
- func (v Vec2) DistSquared(v1 Vec2) float64
- func (v Vec2) Dot(v1 Vec2) float64
- func (v Vec2) Len() int
- func (v Vec2) Max(v1 Vec2) Vec2
- func (v Vec2) Min(v1 Vec2) Vec2
- func (v Vec2) Norm() float64
- func (v Vec2) Normalize() Vec2
- func (v Vec2) ProjectOut(v1 Vec2) Vec2
- func (v Vec2) Scale(f float64) Vec2
- func (v Vec2) Sub(v1 Vec2) Vec2
- func (v Vec2) Sum() float64
- func (v Vec2) WithDim(i int, x float64) Vec2
- func (v Vec2) Zeros() Vec2
- type Vec3
- func (v Vec3) Add(v1 Vec3) Vec3
- func (v Vec3) At(i int) float64
- func (v Vec3) Cross(v1 Vec3) Vec3
- func (v Vec3) Dist(v1 Vec3) float64
- func (v Vec3) DistSquared(v1 Vec3) float64
- func (v Vec3) Dot(v1 Vec3) float64
- func (v Vec3) Len() int
- func (v Vec3) Max(v1 Vec3) Vec3
- func (v Vec3) Min(v1 Vec3) Vec3
- func (v Vec3) Norm() float64
- func (v Vec3) Normalize() Vec3
- func (v Vec3) OrthoBasis() (Vec3, Vec3)
- func (v Vec3) ProjectOut(v1 Vec3) Vec3
- func (v Vec3) Scale(f float64) Vec3
- func (v Vec3) Sub(v1 Vec3) Vec3
- func (v Vec3) Sum() float64
- func (v Vec3) WithDim(i int, x float64) Vec3
- func (v Vec3) Zeros() Vec3
- type Vec4
- func (v Vec4) Add(v1 Vec4) Vec4
- func (v Vec4) At(i int) float64
- func (v Vec4) Dist(v1 Vec4) float64
- func (v Vec4) DistSquared(v1 Vec4) float64
- func (v Vec4) Dot(v1 Vec4) float64
- func (v Vec4) Len() int
- func (v Vec4) Max(v1 Vec4) Vec4
- func (v Vec4) Min(v1 Vec4) Vec4
- func (v Vec4) Norm() float64
- func (v Vec4) Normalize() Vec4
- func (v Vec4) OrthoBasis() (Vec4, Vec4, Vec4)
- func (v Vec4) ProjectOut(v1 Vec4) Vec4
- func (v Vec4) Scale(f float64) Vec4
- func (v Vec4) Sub(v1 Vec4) Vec4
- func (v Vec4) Sum() float64
- func (v Vec4) WithDim(i int, x float64) Vec4
- func (v Vec4) Zeros() Vec4
- type Vector
Constants ¶
const DefaultGSSIters = 64
Variables ¶
This section is empty.
Functions ¶
Types ¶
type BiCGSTAB ¶ added in v0.4.0
BiCGSTAB implements the biconjugate gradient stabilized method for inverting a specific large asymmetric matrix.
This can be instantiated with NewBiCGSTAB() and then iterated via the Iter() method, which returns the current solution.
func NewBiCGSTAB ¶ added in v0.4.0
NewBiCGSTAB initializes the BiCGSTAB solver for the given linear system, as implement by op.
If initGuess is non-nil, it is the initial solution.
type BiCGSTABSolver ¶ added in v0.4.0
type BiCGSTABSolver struct { MaxIters int // If either MSE or MAE go below these values, the // optimization will terminate early. MSETolerance float64 MAETolerance float64 }
BiCGSTABSolver implements LargeLinearSolver using the BiCGSTAB algorithm with a specified stopping condition.
The solver is configured using a tolerance on either MSE, MAE, iterations, or a combination thereof. At least one stopping criterion must be provided.
func (*BiCGSTABSolver) SolveLinearSystem ¶ added in v0.4.0
func (b *BiCGSTABSolver) SolveLinearSystem(op func(v Vec) Vec, x, initGuess Vec) Vec
SolveLinearSystem iteratively runs BiCGSTAB until a stopping criterion is met.
type FiniteVector ¶ added in v0.4.2
type FiniteVector[T any] interface { Vector[T] Len() int WithDim(i int, value float64) T At(i int) float64 }
A FiniteVector is a Vector with the extra constraint that it has a finite number of dimensions that can be accessed independently.
type GridSearch2D ¶ added in v0.4.0
type GridSearch2D struct { // Number of points to test along x and y axes. XStops int YStops int // Number of times to recursively search around the // best point in the grid. Recursions int }
A GridSearch2D implements 2D grid search for minimizing or maximizing 2D functions.
type GridSearch3D ¶ added in v0.4.2
type GridSearch3D struct { // Number of points to test along x, y, and z axes. XStops int YStops int ZStops int // Number of times to recursively search around the // best point in the grid. Recursions int }
A GridSearch3D is like GridSearch2D, but for searching over a 3D volume instead of a 2D area.
type KMeans ¶ added in v0.3.3
type KMeans[T Vector[T]] struct { // Centers is the K-means cluster centers. Centers []T // Data stores all of the data points. Data []T }
KMeans stores the (intermediate) results of K-means clustering.
This type can be created using NewKMeans() to create an initial clustering using K-means++. Then, Iterate() can be called repeatedly to refine the clustering as desired.
type LargeLinearSolver ¶ added in v0.4.0
type LargeLinearSolver interface { // SolveLinearSystem applies the inverse of a linear // operator to the vector b. // // The initGuess argument can be ignored, but may be // provided as a starting point for the solver. SolveLinearSystem(op func(v Vec) Vec, b, initGuess Vec) Vec }
A LargeLinearSolver provides numerical solutions to linear systems that are implemented as functions on vectors.
type LineSearch ¶ added in v0.4.0
type LineSearch struct { // Number of points to test along the input space. Stops int // Number of times to recursively search around the // best point on the line. Recursions int }
A LineSearch implements a 1D line search for minimizing or maximizing 1D functions.
The LineSearch is performed by first sampling coarsely along the function, and then iteratively sampling more densely around the current optimal solution.
type Matrix2 ¶ added in v0.3.4
type Matrix2 [4]float64
Matrix2 is a 2x2 matrix, stored in row-major order.
func NewMatrix2Columns ¶ added in v0.3.4
NewMatrix2Columns creates a Matrix2 with the given coordinates as column entries.
func NewMatrix2Rotation ¶ added in v0.3.4
NewMatrix2Rotation creates a rotation matrix that rotates column vectors by theta.
func (*Matrix2) Eigenvalues ¶ added in v0.3.4
func (m *Matrix2) Eigenvalues() [2]complex128
Eigenvalues computes the eigenvalues of the matrix.
There may be a repeated eigenvalue, but for numerical reasons two are always returned.
func (*Matrix2) InvertInPlace ¶ added in v0.3.4
func (m *Matrix2) InvertInPlace()
InvertInPlace moves the inverse of m into m without causing any new allocations.
func (*Matrix2) InvertInPlaceDet ¶ added in v0.3.4
InvertInPlaceDet is an optimization for InvertInPlace when the determinant has been pre-computed.
func (*Matrix2) MulColumn ¶ added in v0.3.4
MulColumn multiplies the matrix m by a column vector represented by c.
func (*Matrix2) MulColumnInv ¶ added in v0.3.4
MulColumnInv multiplies the inverse of m by the column c, given the determinant of m.
func (*Matrix2) SVD ¶ added in v0.3.4
SVD computes the singular value decomposition of the matrix.
It populates matrices u, s, and v, such that
m = u*s*v.Transpose()
The singular values in s are sorted largest to smallest.
type Matrix3 ¶ added in v0.3.4
type Matrix3 [9]float64
Matrix3 is a 3x3 matrix, stored in row-major order.
func NewMatrix3Columns ¶ added in v0.3.4
NewMatrix3Columns creates a Matrix3 with the given coordinates as column entries.
func NewMatrix3Rotation ¶ added in v0.3.4
NewMatrix3Rotation creates a 3D rotation matrix. Points are rotated around the given vector in a right-handed direction.
The axis is assumed to be normalized. The angle is measured in radians.
func (*Matrix3) Eigenvalues ¶ added in v0.3.4
func (m *Matrix3) Eigenvalues() [3]complex128
Eigenvalues computes the eigenvalues of the matrix.
There may be repeated eigenvalues, but for numerical reasons three are always returned.
func (*Matrix3) InvertInPlace ¶ added in v0.3.4
func (m *Matrix3) InvertInPlace()
InvertInPlace moves the inverse of m into m without causing any new allocations.
func (*Matrix3) InvertInPlaceDet ¶ added in v0.3.4
InvertInPlaceDet is an optimization for InvertInPlace when the determinant has been pre-computed.
func (*Matrix3) MulColumn ¶ added in v0.3.4
MulColumn multiplies the matrix m by a column vector represented by c.
func (*Matrix3) MulColumnInv ¶ added in v0.3.4
MulColumnInv multiplies the inverse of m by the column c, given the determinant of m.
func (*Matrix3) SVD ¶ added in v0.3.4
SVD computes the singular value decomposition of the matrix.
It populates matrices u, s, and v, such that
m = u*s*v.Transpose()
The singular values in s are sorted largest to smallest.
type Matrix4 ¶ added in v0.3.4
type Matrix4 [16]float64
func NewMatrix4Columns ¶ added in v0.3.4
NewMatrix4Columns creates a Matrix4 from the columns.
func NewMatrix4Identity ¶ added in v0.3.4
func NewMatrix4Identity() *Matrix4
NewMatrix4Identity creates an identity matrix.
func (*Matrix4) CharPoly ¶ added in v0.3.4
func (m *Matrix4) CharPoly() Polynomial
CharPoly computes the characteristic polynomial of the matrix.
func (*Matrix4) SVD ¶ added in v0.3.4
SVD computes the singular value decomposition of the matrix.
It populates matrices u, s, and v, such that
m = u*s*v.Transpose()
The singular values in s are sorted largest to smallest.
type Polynomial ¶
type Polynomial []float64
A Polynomial is an equation of the form
a0 + a1*x + a2*x^2 + a3*x^3 + ...
Here, the polynomial is represented as an array of [a0, a1, ...].
func (Polynomial) Add ¶
func (p Polynomial) Add(p1 Polynomial) Polynomial
Add returns the sum of p and p1.
func (Polynomial) Derivative ¶
func (p Polynomial) Derivative() Polynomial
Derivative computes the derivative of the polynomial.
func (Polynomial) Eval ¶
func (p Polynomial) Eval(x float64) float64
Eval evaluates the polynomial at the given x value.
func (Polynomial) IterRealRoots ¶
func (p Polynomial) IterRealRoots(f func(x float64) bool)
IterRealRoots iterates over the real roots of p. This is similar to RealRoots(), but allows the caller to stop iteration early be returning false.
func (Polynomial) Mul ¶
func (p Polynomial) Mul(p1 Polynomial) Polynomial
Mul computes the product of two polynomials.
func (Polynomial) RealRoots ¶
func (p Polynomial) RealRoots() []float64
RealRoots computes the real roots of p, i.e. values of X such that p(x) = 0. The result may have duplicates since roots can be repeated.
If the polynomial has an infinite number of roots, one NaN root is returned.
func (Polynomial) Scale ¶ added in v0.3.0
func (p Polynomial) Scale(c float64) Polynomial
Scale returns a new polynomial with the coefficients of p multiplied by c.
func (Polynomial) String ¶
func (p Polynomial) String() string
String returns a string representation of p.
type RecursiveLineSearch ¶ added in v0.4.2
type RecursiveLineSearch[T FiniteVector[T]] struct { LineSearch LineSearch }
A RecursiveLineSearch searches for an optimum in an N-dimensional space by recursively applying line searches along each axis.
func (*RecursiveLineSearch[T]) Maximize ¶ added in v0.4.2
func (r *RecursiveLineSearch[T]) Maximize(min, max T, f func(T) float64) (solution T, fVal float64)
func (*RecursiveLineSearch[T]) Minimize ¶ added in v0.4.2
func (r *RecursiveLineSearch[T]) Minimize(min, max T, f func(T) float64) (solution T, fVal float64)
type SparseCholesky ¶
type SparseCholesky struct {
// contains filtered or unexported fields
}
SparseCholesky is a sparse LU decomposition of a symmetric matrix.
Once instantiated, this object can be used to quickly apply the inverse of a matrix to vectors.
func NewSparseCholesky ¶
func NewSparseCholesky(mat *SparseMatrix) *SparseCholesky
func (*SparseCholesky) ApplyInverseVec2 ¶ added in v0.4.0
func (s *SparseCholesky) ApplyInverseVec2(x []Vec2) []Vec2
ApplyInverseVec2 computes (A^-1*x, A^-1*y).
func (*SparseCholesky) ApplyInverseVec3 ¶
func (s *SparseCholesky) ApplyInverseVec3(x []Vec3) []Vec3
ApplyInverseVec3 computes (A^-1*x, A^-1*y, A^-1*z).
func (*SparseCholesky) ApplyVec2 ¶ added in v0.4.0
func (s *SparseCholesky) ApplyVec2(x []Vec2) []Vec2
ApplyVec2 computes (A*x, A*y).
func (*SparseCholesky) ApplyVec3 ¶
func (s *SparseCholesky) ApplyVec3(x []Vec3) []Vec3
ApplyVec3 computes (A*x, A*y, A*z).
type SparseMatrix ¶
type SparseMatrix struct {
// contains filtered or unexported fields
}
A SparseMatrix is a square matrix where entries can be set to non-zero values in any order, and not all entries must be set.
func NewSparseMatrix ¶
func NewSparseMatrix(size int) *SparseMatrix
func (*SparseMatrix) Apply ¶ added in v0.4.0
func (s *SparseMatrix) Apply(x Vec) Vec
Apply computes A*x.
func (*SparseMatrix) ApplyVec2 ¶ added in v0.4.0
func (s *SparseMatrix) ApplyVec2(x []Vec2) []Vec2
ApplyVec2 computes (A*x, A*y).
func (*SparseMatrix) ApplyVec3 ¶
func (s *SparseMatrix) ApplyVec3(x []Vec3) []Vec3
ApplyVec3 computes (A*x, A*y, A*z).
func (*SparseMatrix) Iterate ¶
func (s *SparseMatrix) Iterate(row int, f func(col int, x float64))
Iterate loops through the non-zero entries in a row.
func (*SparseMatrix) Permute ¶
func (s *SparseMatrix) Permute(perm []int) *SparseMatrix
Permute permutes the rows and columns by perm, where perm is the result of applying the permutation to the list [0...n-1].
func (*SparseMatrix) RCM ¶
func (s *SparseMatrix) RCM() []int
RCM computes the reverse Cuthill-McKee permutation for the matrix.
func (*SparseMatrix) Set ¶
func (s *SparseMatrix) Set(row, col int, x float64)
Set adds an entry to the matrix.
The entry should not already be set. If it is, this will result in undefined behavior.
func (*SparseMatrix) Transpose ¶ added in v0.4.0
func (s *SparseMatrix) Transpose() *SparseMatrix
Transpose computes the matrix transpose of s.
type Vec ¶
type Vec []float64
Vec is a vector of arbitrary dimension.
func (Vec) DistSquared ¶ added in v0.3.3
DistSquared computes ||v-v1||^2.
func (Vec) ProjectOut ¶
ProjectOut projects the direction v1 out of v.
type Vec2 ¶ added in v0.3.4
type Vec2 [2]float64
A Vec2 is a 2-dimensional tuple of floats.
func NewVec2RandomNormal ¶ added in v0.3.4
func NewVec2RandomNormal() Vec2
NewVec2RandomNormal creates a random normal vector.
func NewVec2RandomUnit ¶ added in v0.3.4
func NewVec2RandomUnit() Vec2
NewVec2RandomUnit creates a unit-length Vec2 in a random direction.
func (Vec2) DistSquared ¶ added in v0.4.0
DistSquared computes ||v-v1||^2.
func (Vec2) ProjectOut ¶ added in v0.3.4
ProjectOut projects the v1 direction out of v.
type Vec3 ¶
type Vec3 [3]float64
A Vec3 is a 3-dimensional tuple of floats.
func LeastSquares3 ¶ added in v0.3.4
LeastSquares3 computes the least squares solution to the equation Ax = b, where A is a matrix with rows vs, and b is a column matrix.
The epsilon argument is a lower bound for singular values in the psuedoinverse.
func LeastSquaresReg3 ¶ added in v0.3.4
LeastSquaresReg3 is like LeastSquares3, but uses ridge regression with a penalty equal to lambda.
func NewVec3RandomNormal ¶
func NewVec3RandomNormal() Vec3
NewVec3RandomNormal creates a random normal vector.
func NewVec3RandomUnit ¶ added in v0.3.4
func NewVec3RandomUnit() Vec3
NewVec3RandomUnit creates a unit-length Vec3 in a random direction.
func (Vec3) DistSquared ¶ added in v0.4.0
DistSquared computes ||v-v1||^2.
func (Vec3) OrthoBasis ¶ added in v0.3.4
OrthoBasis creates two unit vectors which are orthogonal to v and to each other.
If v is axis-aligned, the other vectors will be as well.
func (Vec3) ProjectOut ¶ added in v0.3.4
ProjectOut projects the v1 direction out of v.
type Vec4 ¶ added in v0.3.4
type Vec4 [4]float64
A Vec4 is a 4-dimensional tuple of floats.
func NewVec4RandomNormal ¶ added in v0.3.4
func NewVec4RandomNormal() Vec4
NewVec4RandomNormal creates a random normal vector.
func NewVec4RandomUnit ¶ added in v0.3.4
func NewVec4RandomUnit() Vec4
NewVec4RandomUnit creates a unit-length Vec3 in a random direction.
func (Vec4) DistSquared ¶ added in v0.4.0
DistSquared computes ||v-v1||^2.
func (Vec4) OrthoBasis ¶ added in v0.3.4
OrthoBasis creates three unit vectors which are orthogonal to v and to each other.
If v is axis-aligned, the other vectors will be as well.
The behavior of this method is undefined if v is zero.
func (Vec4) ProjectOut ¶ added in v0.3.4
ProjectOut projects the v1 direction out of v.
type Vector ¶ added in v0.4.0
type Vector[T any] interface { Zeros() T Add(T) T Sub(T) T Scale(float64) T DistSquared(T) float64 Dist(T) float64 Norm() float64 Min(T) T Max(T) T }
The Vector interface type is meant to be used in type constraints as `T Vector[T]`. In this way, Vec2, Vec3, Vec4, and Vec all implement this interface.