class LA::GeneralMatrix(T)

Overview

generic matrix, heap-allocated

Data are stored in column-major as this is a storage used by LAPACK

See SUPPORTED_TYPES for supported types

Included Modules

Defined in:

linalg/cholesky.cr
linalg/eig.cr
linalg/expm.cr
linalg/linalg.cr
linalg/lu.cr
linalg/mult.cr
linalg/qr.cr
linalg/rq_lq_ql.cr
linalg/schur.cr
matrix/general_matrix.cr

Constructors

Class Method Summary

Instance Method Summary

Instance methods inherited from class LA::Matrix(T)

*(k : Number)
*(k : Complex)
*(m : Matrix(T))
*
, **(other : Int) **, +(k : Number)
+(k : Complex)
+(m : Matrix(T))
+
, -(k : Number | Complex)
-(m : Matrix(T))
-
-
, /(k : Number | Complex) /, ==(other) ==, [](i : Int32, j : Int32)
[](arows : Range(Int32 | Nil, Int32 | Nil), acolumns : Range(Int32 | Nil, Int32 | Nil))
[](row : Int32, acolumns : Range(Int32 | Nil, Int32 | Nil))
[](arows : Range(Int32 | Nil, Int32 | Nil), column : Int32)
[]
, []=(i : Int32, j : Int32, value)
[]=(arows : Range(Int32, Int32), acolumns : Range(Int32, Int32), value)
[]=(row : Int32, acolumns : Range(Int32, Int32), value)
[]=(nrows : Range(Int32, Int32), column : Int32, value)
[]=
, add(m : Matrix, *, alpha = 1, beta = 1) add, add!(k : Number, m : Matrix) add!, almost_eq(other : Matrix(T), eps)
almost_eq(other : Matrix(T))
almost_eq
, assume!(flag : MatrixFlags, value : Bool = true) assume!, balance(*, permute = true, scale = true, separate = false) balance, cat(other : Matrix(T), axis : Axis) cat, cholesky(*, lower = false, dont_clean = false) cholesky, chop(eps = self.tolerance) chop, clear_flags clear_flags, columns columns, conjt conjt, conjt! conjt!, conjtranspose conjtranspose, coshm coshm, cosm cosm, det det, detect(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance) detect, detect?(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance) detect?, diag(offset = 0) diag, each(*, all = false, &) each, each_index(*, all = false, &) each_index, each_lower(*, diagonal = true, all = false, &) each_lower, each_upper(*, diagonal = true, all = false, &) each_upper, each_with_index(*, all = false, &) each_with_index, eigs(*, left = false)
eigs(*, need_left : Bool, need_right : Bool)
eigs(*, b : self, need_left = false, need_right = false)
eigs
, eigvals eigvals, expm expm, flags : MatrixFlags flags, hcat(other) hcat, hessenberg
hessenberg(*, calc_q = false)
hessenberg
, inspect(io) inspect, inv inv, kron(b : Matrix(T)) kron, lq lq, lq_r lq_r, lu lu, lu_factor : LUMatrix(T) lu_factor, map(&block : T -> U) forall U map, map!(&) map!, map_with_index(&block : T -> U) forall U map_with_index, map_with_index!(&) map_with_index!, max(axis : Axis) max, min(axis : Axis) min, ncolumns : Int32 ncolumns, norm(kind : MatrixNorm = MatrixNorm::Frobenius) norm, nrows : Int32 nrows, pinv pinv, product(axis : Axis) product, ql ql, ql_r ql_r, qr(*, pivoting = false) qr, qr_r(*, pivoting = false) qr_r, qr_raw(*, pivoting = false) qr_raw, qz(b : self) qz, reduce(axis : Axis, initial, &) reduce, repmat(arows, acolumns) repmat, rows rows, rq rq, rq_r rq_r, save_csv(filename) save_csv, scale!(k : Number | Complex) scale!, schur schur, shape shape, shape_str shape_str, sinhm sinhm, sinm sinm, size : Tuple(Int32, Int32) size, solve(b : self)
solve(b : GeneralMatrix(T), *, overwrite_b = false)
solve
, square? square?, sum(axis : Axis) sum, svd svd, svdvals svdvals, t t, t! t!, tanhm tanhm, tanm tanm, to_custom(io, prefix, columns_separator, rows_separator, postfix)
to_custom(prefix, columns_separator, rows_separator, postfix)
to_custom
, to_general to_general, to_imag to_imag, to_matlab(io)
to_matlab
to_matlab
, to_real to_real, to_s(io) to_s, to_unsafe to_unsafe, tolerance tolerance, trace trace, transpose transpose, tril(k = 0) tril, tril!(k = 0) tril!, triu(k = 0) triu, triu!(k = 0) triu!, vcat(other) vcat

Class methods inherited from class LA::Matrix(T)

arange(start_val : T, end_val : T, delta = 1.0) arange, block_diag(*args) block_diag, circulant(c) circulant, column(values) column, companion(a) companion, dft(n, scale : DFTScale = DFTScale::None) dft, diag(nrows, ncolumns, value : Number | Complex)
diag(nrows, ncolumns, values)
diag(values)
diag(nrows, ncolumns, &)
diag
, eye(n) eye, fiedler(values) fiedler, from_custom(str : String, prefix, columns_separator, rows_separator, postfix)
from_custom(io, prefix, columns_separator, rows_separator, postfix)
from_custom
, from_matlab(s) from_matlab, hadamard(n) hadamard, hankel(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil) hankel, helmert(n, full = false) helmert, hilbert(n) hilbert, identity(n) identity, invhilbert(n) invhilbert, invpascal(n, kind : PascalKind = PascalKind::Symmetric) invpascal, kron(a, b) kron, leslie(f, s) leslie, load_csv(filename) load_csv, ones(nrows, ncolumns) ones, pascal(n, kind : PascalKind = PascalKind::Symmetric) pascal, rand(nrows, ncolumns, rng : Random | Nil = nil) rand, repmat(a : Matrix(T), nrows, ncolumns) repmat, row(values) row, toeplitz(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil) toeplitz, tri(nrows, ncolumns, k = 0) tri, zeros(nrows, ncolumns) zeros

Macros inherited from class LA::Matrix(T)

lapack(name, *args, worksize = nil) lapack, lapack_util(name, worksize, *args) lapack_util

Instance methods inherited from module Enumerable(T)

product(initial : Complex)
product(initial : Complex, &)
product

Constructor Detail

def self.new(nrows : Int32, ncolumns : Int32, values : Indexable, col_major = false, flags : LA::MatrixFlags = MatrixFlags::None) #

Creates matrix with given size and populate elements from values

if col_major is true, values content is just copied to #raw, otherwise a conversion from row-major form is performed Example:

values = [1, 2, 3, 4]
a = GMat.new(2, 2, values)
a.to_aa # => [[1,2],[3,4]]
b = GMat.new(2, 2, values, col_major: true)
b.to_aa # => [[1,3],[2,4]]

[View source]
def self.new(nrows : Int32, ncolumns : Int32, flags : LA::MatrixFlags = MatrixFlags::None) #

Creates zero-initialized matrix of given size

Example: LA::GMat.new(4,4)


[View source]
def self.new(nrows : Int32, ncolumns : Int32, flags : LA::MatrixFlags = MatrixFlags::None, &) #

Creates matrix of given size and then call block to initialize each element

Example: LA::GMat.new(4,4){|i,j| i+j }


[View source]
def self.new(values : Indexable) #

Creates matrix from any Indexable of Indexables

Example:

m = GMat.new([
  [1, 2, 3],
  [4, 5, 6],
  [7, 8, 9],
  [10, 11, 12],
])

[View source]
def self.new(matrix : Matrix) #

Creates matrix with same content as another matrix


[View source]

Class Method Detail

def self.[](*args) #

Alias for #new


[View source]
def self.columns(*args) #

Creates matrix from a number of columns

Example:

a = GMat.columns([1, 2, 3, 4], [5, 6, 7, 8])
a.to_aa # => [[1, 5], [2, 6], [3, 7], [4, 8]]

[View source]
def self.diag(nrows, ncolumns, values) #

Returns diagonal matrix of given size with diagonal elements taken from array values


[View source]
def self.diag(nrows, ncolumns, &) #

Returns diagonal matrix of given size with diagonal elements equal to block value


[View source]
def self.rows(*args) #

Creates matrix from a number of rows

Example:

a = GMat.rows([1, 2, 3, 4], [5, 6, 7, 8])
a.to_aa # => [[1,2,3,4], [5,6,7,8]]

[View source]

Instance Method Detail

def ==(other : self) #

[View source]
def abs(kind : MatrixNorm = MatrixNorm::Frobenius) #

Alias for #norm.


[View source]
def add!(m : Matrix, *, alpha = 1, beta = 1) #

Perform inplace linear combination with matrix m a.add!(m, alpha: alpha, beta: beta) is equal to a = alpha * a + beta * m, but faster as no new matrix is allocated


[View source]
def add_mult(a, b : Matrix(T), *, alpha = 1.0, beta = 1.0) #

performs c = alphaab + beta*c (BLAS routines gemm/symm/hemm)


[View source]
def balance!(*, permute = true, scale = true, separate = false) #

Balances a square matrix in-place to improve eigenvalue accuracy.

Arguments:

  • permute (Bool) : If true, performs permutations. Default: true.
  • scale (Bool) : If true, performs scaling. Default: true.
  • separate (Bool) : If true, returns scaling factors separately. Default: false.

Returns:

  • GeneralMatrix(T) : If separate is true, returns scaling factors; otherwise, returns diagonal matrix.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routine:

  • Uses gebal (balance matrix).

[View source]
def cat!(other : Matrix(T), dimension) #

Concatenate matrix adding another matrix by dimension Axis::Rows (horizontal) or Axis::Columns (vertical)


[View source]
def cho_solve(b : self, *, overwrite_b = false) #

Solves A * X = B given Cholesky factorization of A (A = L*L' or U'*U).

Arguments:

  • b (GeneralMatrix(T)) : Right-hand side matrix.
  • overwrite_b (Bool) : Allows overwriting b with solution. Default: false.

Returns:

  • GeneralMatrix(T) : Solution matrix X.

Raises:

  • ArgumentError : If #nrows mismatch, self is not square, or not triangular.

LAPACK Routine:

  • Uses potrs (positive-definite solve).

[View source]
def cholesky!(*, lower = false, dont_clean = false) #

Performs in-place Cholesky decomposition.

Arguments:

  • lower (Bool) : If true, computes lower triangular factor (L). Otherwise, upper (U). Default: false.
  • dont_clean (Bool) : If true, leaves the unused triangle unmodified. Default: false.

Returns:

  • self : On successful decomposition, self becomes the triangular factor.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routine:

  • Uses potrf (real/complex positive-definite factorization).

[View source]
def clone #

returns copy of matrix


[View source]
def conjtranspose! #

Conjurgate transposes matrix inplace

Currently, transpose of non-square matrix still allocates temporary buffer


[View source]
def det(*, overwrite_a = false) #

Calculates the determinant of a square matrix.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • T : Determinant value.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routine:

  • Uses getrf (LU factorization).

[View source]
def dup #

returns copy of matrix


[View source]
def eigs(*, left = false, overwrite_a = false) #

[View source]
def eigs(*, need_left : Bool, need_right : Bool, overwrite_a = false) #

Computes eigenvalues and optionally left and right eigenvectors.

Arguments:

  • need_left: Whether to compute left eigenvectors.
  • need_right: Whether to compute right eigenvectors.
  • overwrite_a: If true, allows overwriting matrix self.

Returns: A tuple {values, left_vectors, right_vectors} where:

  • values is an Array(T) of eigenvalues (or Array(Complex) if complex eigenvalues arise in real matrices).
  • left_vectors is a GeneralMatrix(T)? containing left eigenvectors if need_left == true, else nil.
  • right_vectors is a GeneralMatrix(T)? containing right eigenvectors if need_right == true, else nil.

For Hermitian (complex) or symmetric (real) matrices, eigenvalues are real, and eigenvectors are orthonormal.

Exceptions:

  • ArgumentError: if the matrix is not square.

Called LAPACK routines:

  • heevr — for complex Hermitian matrices.
  • syevr — for real symmetric matrices.
  • geev — for general complex or real non-symmetric/Hermitian matrices.

[View source]
def eigs(*, b : GeneralMatrix(T), need_left : Bool, need_right : Bool, overwrite_a = false, overwrite_b = false) #

Computes generalized eigenvalues and optionally left and right eigenvectors.

Arguments:

  • b: Right-hand matrix in the generalized eigenvalue problem.
  • need_left: Whether to compute left eigenvectors.
  • need_right: Whether to compute right eigenvectors.
  • overwrite_a: If true, allows overwriting matrix self.
  • overwrite_b: If true, allows overwriting matrix b.

Returns: A tuple {alpha, beta, left_vectors, right_vectors} where:

  • alpha is an Array(T) of generalized eigenvalue numerators (or Array(Complex) if complex values arise).
  • beta is an Array(T) of generalized eigenvalue denominators (λ = alpha[i] / beta[i]).
  • left_vectors is a GeneralMatrix(T)? containing left eigenvectors if need_left == true, else nil.
  • right_vectors is a GeneralMatrix(T)? containing right eigenvectors if need_right == true, else nil.

For Hermitian/symmetric A and positive definite B, the problem is solved using a specialized algorithm ensuring real eigenvalues.

Exceptions:

  • ArgumentError: if self is not square or if b does not have the same size as self.

Called LAPACK routines:

  • hegvd — for complex Hermitian-definite generalized problems.
  • sygvd — for real symmetric-definite generalized problems.
  • ggev — for general complex or real non-symmetric/Hermitian problems.

[View source]
def eigvals(*, overwrite_a = false) #

[View source]
def expm(*, schur_fact = false) #

Computes matrix exponential

It exploits triangularity (if any) of A. if schur_fact is true, function uses an initial transformation to Schur form.

Reference: A. H. Al-Mohy and N. J. Higham, A New Scaling and Squaring Algorithm for the Matrix Exponential, SIAM J. Matrix Anal. Appl. 31(3): 970-989, 2009. Awad H. Al-Mohy and Nicholas J. Higham, April 20, 2010.


[View source]
def flags : MatrixFlags #

Matrix flags (see MatrixFlags for description)


[View source]
def flags=(flags : MatrixFlags) #

Matrix flags (see MatrixFlags for description)


[View source]
def hcat!(other) #

Concatenate matrix adding another matrix horizontally (so they form a row)


[View source]
def hessenberg! #

Computes the Hessenberg form in-place (without Q).

Returns:

  • self : Hessenberg matrix.

[View source]
def hessenberg!(*, calc_q = false) #

Computes the Hessenberg form of a matrix in-place.

Arguments:

  • calc_q (Bool) : If true, also computes the orthogonal matrix Q. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T)) : Hessenberg matrix and optionally Q.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routines:

  • gebal (balance)
  • gehrd (Hessenberg reduction)
  • orghr (generate Q)

[View source]
def inv! #

Calculate matrix inversion in-place.

Method selects optimal algorithm depending on MatrixFlags:

  • #transpose! returned if matrix is orthogonal.
  • trtri used if matrix is triangular.
  • potrf + potri used if matrix is positive definite.
  • hetrf + hetri used if matrix is hermitian (complex only).
  • sytrf + sytri used if matrix is symmetric.
  • getrf + getri used otherwise.

Returns:

  • self : On successful inversion, self becomes the inverse.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routines:

  • trtri (triangular inverse)
  • potrf/potri (positive definite)
  • hetrf/hetri (hermitian)
  • sytrf/sytri (symmetric)
  • getrf/getri (general)

[View source]
def lq(*, overwrite_a = false) #

Computes LQ decomposition.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T)) : L and Q factors.

LAPACK Routines:

  • gelqf (LQ factorization)
  • orglq (generate Q from elementary reflectors)

[View source]
def lstsq(b : self, method : LSMethod = LSMethod::Auto, *, overwrite_a = false, overwrite_b = false, cond = -1) #

Solves the linear least squares problem with multiple method options.

Arguments:

  • b (GeneralMatrix(T)) : Right-hand side matrix.
  • method (LSMethod) : Algorithm to use. Default: LSMethod::Auto.
  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.
  • overwrite_b (Bool) : If true, allows overwriting b. Default: false.
  • cond (Number) : Cutoff for rank determination. Default: -1 (machine precision).

Returns:

  • Tuple(GeneralMatrix(T), Int32, Array(T)) : Solution X, effective rank, singular values (if applicable).

Raises:

  • ArgumentError : If #nrows mismatch.

LAPACK Routines:

  • gels (QR least squares)
  • gelsd (SVD-based least squares)
  • gelsy (QR with pivoting)

[View source]
def lu(*, overwrite_a = false) #

Compute pivoted LU decomposition of a matrix

If you call p,l,u = a.lu for a matrix m*n a then

  • p will be m*m permutation matrix
  • l will be m*k lower triangular matrix with unit diagonal. K = min(M, N)
  • u will be k*n upper triangular matrix
  • (p*l*u) will be equal to a (within calculation tolerance)

if overwrite_a is true, a will be overriden in process of calculation

Uses getrf LAPACK routine


[View source]
def lu_factor! : LUMatrix(T) #

Compute pivoted LU decomposition of a matrix and returns it in a compact form, useful for solving linear equation systems.

Overrides source matrix in a process of calculation

Uses getrf LAPACK routine


[View source]
def ncolumns : Int32 #

Count of columns in matrix


[View source]
def norm(kind : MatrixNorm = MatrixNorm::Frobenius) #

Computes the norm of the matrix.

Method selects optimal algorithm depending on MatrixFlags:

  • lantr used for triangular matrices.
  • lanhe used for hermitian matrices (complex only).
  • lansy used for symmetric matrices.
  • lange used for general matrices.

Arguments:

Returns:

  • T : Norm value.

LAPACK Routines:

  • lantr (triangular matrix norm)
  • lanhe (hermitian matrix norm)
  • lansy (symmetric matrix norm)
  • lange (general matrix norm)

[View source]
def nrows : Int32 #

Count of rows in matrix


[View source]
def ql(*, overwrite_a = false) #

Computes QL decomposition.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T)) : Q and L factors.

LAPACK Routines:

  • geqlf (QL factorization)
  • orgql (generate Q from elementary reflectors)

[View source]
def qr(*, overwrite_a = false, pivoting = false) #

Computes the complete QR decomposition with orthogonal Q.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.
  • pivoting (Bool) : If true, performs QR with column pivoting. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T), Array(Int32)) : Orthogonal matrix Q (m×m or m×n depending on dimensions) Upper triangular matrix R (n×n) Pivot indices (empty if pivoting = false)

Properties:

  • For an m×n matrix A with m ≥ n: A = Q * R
  • For m < n: A = Q * R where Q is m×m and R is m×n
  • When pivoting is enabled: A * P = Q * R where P is the permutation matrix defined by pivot indices

Side Effects:

  • Marks Q as Orthogonal on successful decomposition

LAPACK Routines:

  • geqrf/geqp3 (QR factorization)
  • orgqr (generate Q from elementary reflectors)

[View source]
def qr_r(*, overwrite_a = false, pivoting = false) #

Performs in-place QR factorization, returning only the R factor.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.
  • pivoting (Bool) : If true, uses QR with column pivoting. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), Array(Int32)) : Upper triangular matrix R and pivot indices.

LAPACK Routines:

  • geqrf (QR without pivoting)
  • geqp3 (QR with column pivoting)

[View source]
def qr_raw(*, overwrite_a = false, pivoting = false) #
Description copied from class LA::Matrix(T)

See GeneralMatrix#qr_raw


[View source]
def qz(b : self, overwrite_a = false, overwrite_b = false) #

Performs in-place QZ (generalized Schur) decomposition.

Arguments:

  • b (GeneralMatrix(T)) : Second matrix in the generalized eigenvalue problem.
  • overwrite_a (Bool) : If true, allows overwriting the first matrix. Default: false.
  • overwrite_b (Bool) : If true, allows overwriting the second matrix. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T), GeneralMatrix(T), GeneralMatrix(T)) : A tuple containing:
    • Generalized Schur form of A.
    • Generalized Schur form of B.
    • Left transformation matrix VSL.
    • Right transformation matrix VSR.

Raises:

  • ArgumentError : If either matrix is not square.

LAPACK Routine:

  • Uses gges (generalized Schur decomposition).

[View source]
def rank(eps = self.tolerance, *, method : RankMethod = RankMethod::SVD, overwrite_a = false) #

Determines the effective rank of the matrix.

Arguments:

  • eps (Number) : Tolerance for singular value comparison. Default: self.tolerance.
  • method (RankMethod) : Method to use. Default: RankMethod::SVD.
  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Int32 : Estimated rank.

LAPACK Routines:

  • SVD method uses gesdd (singular values)
  • QRP method uses QR with pivoting

[View source]
def raw : Slice(T) #

Pointer to a raw data


[View source]
def reshape(anrows, ancolumns, col_major = false) #

Returns a matrix with different nrows and ncolumns but same elements (total number of elements must not change)

if col_major is true, just nrows and ncolumns are changed, data kept the same Otherwise, elements are reordered to emulate row-major storage Example:

a = GMat[[1, 2, 3], [4, 5, 6]]
a.reshape(2, 3).to_aa # => [[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]]
b = GMat[[1, 2, 3], [4, 5, 6]]
this is because in-memory order of values is (1, 4, 2, 5, 3, 6)
b.reshape(2, 3, col_major: true).to_aa # => [[1.0, 5.0], [4.0, 3.0], [2.0, 6.0]]

[View source]
def reshape!(anrows : Int32, ancolumns : Int32, col_major = false) #

Changes nrows and ncolumns of matrix (total number of elements must not change)

if col_major is true, just nrows and ncolumns are changed, data kept the same Otherwise, elements are reordered to emulate row-major storage Example:

a = GMat[[1, 2, 3], [4, 5, 6]]
a.reshape!(2, 3)
a.to_aa # => [[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]]
b = GMat[[1, 2, 3], [4, 5, 6]]
b.reshape!(2, 3, col_major: true)
this is because in-memory order of values is (1, 4, 2, 5, 3, 6)
b.to_aa # => [[1.0, 5.0], [4.0, 3.0], [2.0, 6.0]]

[View source]
def resize!(anrows : Int32, ancolumns : Int32) #

Change number of rows and columns in matrix.

if new number is higher zero elements are added, if new number is lower, exceeding elements are lost


[View source]
def rq(*, overwrite_a = false) #

Computes RQ decomposition.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T)) : Q and R factors.

LAPACK Routines:

  • gerqf (RQ factorization)
  • orgrq (generate Q from elementary reflectors)

[View source]
def schur(*, overwrite_a = false) #

Performs in-place Schur decomposition.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting the matrix with its Schur form. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), GeneralMatrix(T)) : A tuple containing:
    • Schur form matrix (upper triangular for complex, quasi-triangular for real).
    • Orthogonal/unitary transformation matrix Z.

Raises:

  • ArgumentError : If matrix is not square.

LAPACK Routine:

  • Uses gees (Schur decomposition).

[View source]
def solve(b : self, *, overwrite_a = false, overwrite_b = false) #

Solves the linear system A * X = B.

Method selects optimal algorithm depending on MatrixFlags:

  • trtrs used if matrix is triangular.
  • posv used if matrix is positive definite.
  • hesv used if matrix is hermitian (complex only).
  • sysv used if matrix is symmetric.
  • gesv used otherwise.

Arguments:

  • b (GeneralMatrix(T)) : Right-hand side matrix.
  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.
  • overwrite_b (Bool) : If true, allows overwriting b. Default: false.

Returns:

  • GeneralMatrix(T) : Solution matrix X.

Raises:

  • ArgumentError : If #nrows mismatch or matrix not square.

LAPACK Routines:

  • trtrs (triangular solve)
  • posv (positive definite solve)
  • hesv (hermitian solve)
  • sysv (symmetric solve)
  • gesv (general solve)

[View source]
def solvels(b : self, *, overwrite_a = false, overwrite_b = false, cond = -1) #

Solves the linear least squares problem using QR factorization.

Arguments:

  • b (GeneralMatrix(T)) : Right-hand side matrix.
  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.
  • overwrite_b (Bool) : If true, allows overwriting b. Default: false.
  • cond (Number) : Cutoff for rank determination. Default: -1 (machine precision).

Returns:

  • GeneralMatrix(T) : Solution matrix X (may include residuals).

Raises:

  • ArgumentError : If #nrows mismatch.

LAPACK Routine:

  • Uses gels (QR least squares).

[View source]
def svd(*, overwrite_a = false) #

Computes the singular value decomposition (SVD) of a matrix.

For an m×n matrix A, returns:

  • u : m×m orthogonal matrix
  • s : Array of singular values of length min(m,n)
  • vt : n×n orthogonal matrix

Such that A = u * diag(s) * vt.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Tuple(GeneralMatrix(T), Array(T), GeneralMatrix(T)) : U, singular values, V^T.

LAPACK Routine:

  • Uses gesdd (divide-and-conquer SVD).

[View source]
def svdvals(*, overwrite_a = false) #

Computes only the singular values of a matrix.

Arguments:

  • overwrite_a (Bool) : If true, allows overwriting self. Default: false.

Returns:

  • Array(T) : Singular values.

LAPACK Routine:

  • Uses gesdd (SVD, singular values only).

[View source]
def to_a(col_major = false) #

Converts matrix to plain array of elements if col_major is true, elements are returned as stored inplace, otherwise row-major storage is emulated


[View source]
def to_aa(col_major = false) #

Converts matrix to array of array of elements if col_major is true, elements are returned as stored inplace, otherwise row-major storage is emulated Example:

a = GMat[[1, 2], [3, 4]]
a.to_aa                  # => [[1.0, 2.0],[3.0, 4.0]]
a.to_aa(col_major: true) # => [[1.0, 3.0],[2.0, 4.0]]

[View source]
def to_unsafe #

Returns pointer to raw data, suitable to e.g. use with LAPACK


[View source]
def tr_mult!(a : Matrix(T), *, alpha = 1.0, left = false) #

performs b = alphaab or b = alphaba (BLAS routine trmm)

a must be a square triangular GeneralMatrix(T)

if left is true, alpha*a*b is calculated, otherwise alpha*b*a


[View source]
def transpose! #

transposes matrix inplace

Currently, transpose of non-square matrix still allocates temporary buffer


[View source]
def unsafe_fetch(i, j) #

returns element at row i and column j, without performing any checks


[View source]
def unsafe_set(i, j, value) #

sets element at row i and column j to value, without performing any checks


[View source]
def vcat!(other) #

Concatenate matrix adding another matrix vertically (so they form a column)


[View source]