lubeck

Lubeck - Linear Algebra

Members

Enums

fastmath
anonymousenum fastmath
Undocumented in source.

Functions

choleskyDecomp
auto choleskyDecomp(Slice!(kind, [2], Iterator) a, char uplo)

Computs Cholesky decomposition of symmetric positive definite matrix 'A'. The factorization has the form: \A = U**T * U, if UPLO = Upper, or \A = L * L**T, if UPLO = Lower Where U is an upper triangular matrix and L is lower triangular.

choleskySolve
auto choleskySolve(Slice!(Canonical, [2], IteratorC) c, Slice!(kindB, n, IteratorB) b, char uplo)

Solves a system of linear equations A * X = B with a symmetric matrix A using the Cholesky factorization: \A = U**T * U or \A = L * L**T computed by choleskyDecomp.

cov
Slice!(Contiguous, [2], BlasType!Iterator*) cov(Slice!(kind, [2], Iterator) matrix)

Covariance matrix.

det
auto det(Slice!(kind, [2], Iterator) a)

Matrix determinant.

detSymmetric
auto detSymmetric(Slice!(kind, [2], Iterator) a, char store)

Matrix determinant.

eigSymmetric
auto eigSymmetric(Slice!(kind, [2], Iterator) a, char store)

Eigenvalues and eigenvectors of symmetric matrix.

inv
auto inv(Slice!(kind, [2], Iterator) a)

Calculates the inverse of a matrix.

ldlDecomp
auto ldlDecomp(Slice!(kind, [2], Iterator) a, char uplo)

Computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method. The for of the factorization is: \A = L*D*L**T Where L is product if permutation and unit lower triangular matrices, and D is symmetric and block diagonal with '1 x 1' and '2 x 2' diagonal blocks.

ldlSolve
auto ldlSolve(Slice!(Canonical, [2], IteratorA) a, Slice!(Contiguous, [1], lapackint*) ipiv, Slice!(kindB, n, IteratorB) b, char uplo)

Solves a system of linear equations \A * X = B with symmetric matrix 'A' using the factorization \A = U * D * U**T, or \A = L * D * L**T computed by ldlDecomp.

luDecomp
auto luDecomp(Slice!(kind, [2], Iterator) a)

Computes LU factorization of a general 'M x N' matrix 'A' using partial pivoting with row interchanges. The factorization has the form: \A = P * L * U Where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

luSolve
auto luSolve(Slice!(Canonical, [2], IteratorLU) lut, Slice!(Contiguous, [1], lapackint*) ipiv, Slice!(kindB, n, IteratorB) b, char trans)

Solves a system of linear equations \A * X = B, or \A**T * X = B with a general 'N x N' matrix 'A' using the LU factorization computed by luDecomp.

mldivide
Slice!(Contiguous, [2], BlasType!(IteratorA, IteratorB)*) mldivide(Slice!(kindA, [2], IteratorA) a, Slice!(kindB, [2], IteratorB) b)
Slice!(Contiguous, [1], BlasType!(IteratorA, IteratorB)*) mldivide(Slice!(kindA, [2], IteratorA) a, Slice!(kindB, [1], IteratorB) b)

Solve systems of linear equations AX = B for X. Computes minimum-norm solution to a linear least squares problem if A is not a square matrix.

mtimes
Slice!(Contiguous, [2], BlasType!(IteratorA, IteratorB)*) mtimes(Slice!(kindA, [2], IteratorA) a, Slice!(kindB, [2], IteratorB) b)

General matrix-matrix multiplication. Allocates result to an uninitialized slice using GC.

mtimes
Slice!(Contiguous, [1], BlasType!(IteratorA, IteratorB)*) mtimes(Slice!(kindA, [2], IteratorA) a, Slice!(kindB, [1], IteratorB) b)

General matrix-matrix multiplication. Allocates result to an uninitialized slice using GC.

mtimes
Slice!(Contiguous, [1], BlasType!(IteratorA, IteratorB)*) mtimes(Slice!(kindB, [1], IteratorB) a, Slice!(kindA, [2], IteratorA) b)

General matrix-matrix multiplication.

mtimes
CommonType!(BlasType!IteratorA, BlasType!IteratorB) mtimes(Slice!(kindB, [1], IteratorB) a, Slice!(kindA, [1], IteratorA) b)

Vector-vector multiplication (dot product).

pca
auto pca(Slice!(kind, [2], Iterator) matrix, Flag!"centerColumns" cc)

Principal component analysis of raw data.

pinv
Slice!(Contiguous, [2], BlasType!Iterator*) pinv(Slice!(kind, [2], Iterator) matrix, BlasType!Iterator tolerance)

Computes Moore-Penrose pseudoinverse of matrix.

qrDecomp
auto qrDecomp(Slice!(kind, [2], Iterator) a)

Computes a QR factorization of matrix 'a'.

qrSolve
auto qrSolve(Slice!(Canonical, [2], IteratorA) a, Slice!(Contiguous, [1], IteratorT) tau, Slice!(kindB, n, IteratorB) b)

Solve the least squares problem: \min ||A * X - B|| Using the QR factorization: \A = Q * R computed by qrDecomp.

svd
auto svd(Slice!(kind, [2], Iterator) matrix, Flag!"slim" slim)

Computes the singular value decomposition.

Structs

EigSymmetricResult
struct EigSymmetricResult(T)

Eigenvalues and eigenvectors POD.

LDLResult
struct LDLResult(T)

Consist LDL factorization;

LUResult
struct LUResult(T)

LUResult consist lu factorization.

PcaResult
struct PcaResult(T)

Principal component analises result.

QRResult
struct QRResult(T)
SvdResult
struct SvdResult(T)
choleskyResult
struct choleskyResult(T)
Undocumented in source.

Meta

Authors

Ilya Yaroshenko, Lars Tandle Kyllingstad (SciD author)