class LA::GeneralMatrix(T)
- LA::GeneralMatrix(T)
- LA::Matrix(T)
- Reference
- Object
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
- LA::DenseMatrix
Defined in:
linalg/cholesky.crlinalg/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
-
.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
-
.new(nrows : Int32, ncolumns : Int32, flags : LA::MatrixFlags = MatrixFlags::None)
Creates zero-initialized matrix of given size
-
.new(nrows : Int32, ncolumns : Int32, flags : LA::MatrixFlags = MatrixFlags::None, &)
Creates matrix of given size and then call block to initialize each element
-
.new(values : Indexable)
Creates matrix from any
Indexable
ofIndexable
s -
.new(matrix : Matrix)
Creates matrix with same content as another matrix
Class Method Summary
-
.[](*args)
Alias for
#new
-
.columns(*args)
Creates matrix from a number of columns
-
.diag(nrows, ncolumns, values)
Returns diagonal matrix of given size with diagonal elements taken from array
values
-
.diag(nrows, ncolumns, &)
Returns diagonal matrix of given size with diagonal elements equal to block value
-
.rows(*args)
Creates matrix from a number of rows
Instance Method Summary
-
#==(other : self)
see
LA::Matrix#==
-
#abs(kind : MatrixNorm = MatrixNorm::Frobenius)
Alias for
#norm
-
#add!(m : Matrix, *, alpha = 1, beta = 1)
Perform inplace linear combination with matrix
m
a.add!(m, alpha: alpha, beta: beta)
is equal toa = alpha * a + beta * m
, but faster as no new matrix is allocated -
#add_mult(a, b : Matrix(T), *, alpha = 1.0, beta = 1.0)
performs c = alphaab + beta*c (BLAS routines
gemm
/symm
/hemm
) - #balance!(*, permute = true, scale = true, separate = false)
-
#cat!(other : Matrix(T), dimension)
Concatenate matrix adding another matrix by
dimension
Axis::Rows
(horizontal) orAxis::Columns
(vertical) - #cho_solve(b : self, *, overwrite_b = false)
- #cholesky!(*, lower = false, dont_clean = false)
-
#clone
returns copy of matrix
-
#conjtranspose!
Conjurgate transposes matrix inplace
-
#det(*, overwrite_a = false)
Calculates determinant for a square matrix
-
#dup
returns copy of matrix
- #eigs(*, left = false, overwrite_a = false)
- #eigs(*, need_left : Bool, need_right : Bool, overwrite_a = false)
-
#eigs(*, b : GeneralMatrix(T), need_left : Bool, need_right : Bool, overwrite_a = false, overwrite_b = false)
generalized eigenvalues problem
- #eigvals(*, overwrite_a = false)
-
#expm(*, schur_fact = false)
Computes matrix exponential
-
#flags : MatrixFlags
Matrix flags (see
MatrixFlags
for description) -
#flags=(flags : MatrixFlags)
Matrix flags (see
MatrixFlags
for description) -
#hcat!(other)
Concatenate matrix adding another matrix horizontally (so they form a row)
- #hessenberg!
- #hessenberg!(*, calc_q = false)
-
#inv!
Calculate matrix inversion inplace Method selects optimal algorithm depending on
MatrixFlags
transpose
returned if matrix is orthogonaltrtri
is used if matrix is triangularpotrf
,potri
are used if matrix is positive definitehetrf
,hetri
are used if matrix is hermitiansytrf
,sytri
are used if matrix is symmetricgetrf
,getri
are used otherwise - #lq(*, overwrite_a = false)
- #lq_r(*, overwrite_a = false)
- #lstsq(b : self, method : LSMethod = LSMethod::Auto, *, overwrite_a = false, overwrite_b = false, cond = -1)
-
#lu(*, overwrite_a = false)
Compute pivoted LU decomposition of a matrix
-
#lu_factor! : LUMatrix(T)
Compute pivoted LU decomposition of a matrix and returns it in a compact form, useful for solving linear equation systems.
-
#ncolumns : Int32
Count of columns in matrix
-
#norm(kind : MatrixNorm = MatrixNorm::Frobenius)
returns matrix norm
-
#nrows : Int32
Count of rows in matrix
- #ql(*, overwrite_a = false)
- #ql_r(*, overwrite_a = false)
- #qr(*, overwrite_a = false, pivoting = false)
- #qr_r(*, overwrite_a = false, pivoting = false)
- #qr_raw(*, overwrite_a = false, pivoting = false)
- #qz(b : self, overwrite_a = false, overwrite_b = false)
-
#rank(eps = self.tolerance, *, method : RankMethod = RankMethod::SVD, overwrite_a = false)
determine effective rank either by SVD method or QR-factorization with pivoting
-
#raw : Slice(T)
Pointer to a raw data
-
#reshape(anrows, ancolumns, col_major = false)
Returns a matrix with different nrows and ncolumns but same elements (total number of elements must not change)
-
#reshape!(anrows : Int32, ancolumns : Int32, col_major = false)
Changes nrows and ncolumns of matrix (total number of elements must not change)
-
#resize!(anrows : Int32, ancolumns : Int32)
Change number of rows and columns in matrix.
- #rq(*, overwrite_a = false)
- #rq_r(*, overwrite_a = false)
- #schur(*, overwrite_a = false)
-
#solve(b : self, *, overwrite_a = false, overwrite_b = false)
Solves matrix equation
self*x = b
and returns x - #solvels(b : self, *, overwrite_a = false, overwrite_b = false, cond = -1)
-
#svd(*, overwrite_a = false)
Calculates singular value decomposition for a matrix
-
#svdvals(*, overwrite_a = false)
Calculates array of singular values for a matrix
-
#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 -
#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]]
-
#to_unsafe
Returns pointer to raw data, suitable to e.g.
-
#tr_mult!(a : Matrix(T), *, alpha = 1.0, left = false)
performs b = alphaab or b = alphaba (BLAS routine
trmm
) -
#transpose!
transposes matrix inplace
-
#unsafe_fetch(i, j)
returns element at row i and column j, without performing any checks
-
#unsafe_set(i, j, value)
sets element at row i and column j to value, without performing any checks
-
#vcat!(other)
Concatenate matrix adding another matrix vertically (so they form a column)
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, 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(&) map, map!(&) map!, map_with_index(&) 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::DEFAULT) 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
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]]
Creates zero-initialized matrix of given size
Example: LA::GMat.new(4,4)
Creates matrix of given size and then call block to initialize each element
Example: LA::GMat.new(4,4){|i,j| i+j }
Creates matrix from any Indexable
of Indexable
s
Example:
m = GMat.new([
[1, 2, 3],
[4, 5, 6],
[7, 8, 9],
[10, 11, 12],
])
Class Method Detail
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]]
Returns diagonal matrix of given size with diagonal elements taken from array values
Returns diagonal matrix of given size with diagonal elements equal to block value
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]]
Instance Method Detail
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
performs c = alphaab + beta*c (BLAS routines gemm
/symm
/hemm
)
Concatenate matrix adding another matrix by dimension
Axis::Rows
(horizontal) or Axis::Columns
(vertical)
Conjurgate transposes matrix inplace
Currently, transpose of non-square matrix still allocates temporary buffer
Calculates determinant for a square matrix
if overwrite_a is true, a will be overriden in process of calculation
Uses getrf
LAPACK routine
generalized eigenvalues problem
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.
Calculate matrix inversion inplace
Method selects optimal algorithm depending on MatrixFlags
transpose
returned if matrix is orthogonal
trtri
is used if matrix is triangular
potrf
, potri
are used if matrix is positive definite
hetrf
, hetri
are used if matrix is hermitian
sytrf
, sytri
are used if matrix is symmetric
getrf
, getri
are used otherwise
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 matrixl
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 toa
(within calculation tolerance)
if overwrite_a is true, a will be overriden in process of calculation
Uses getrf
LAPACK routine
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
returns matrix norm
kind
defines type of norm
Uses LAPACK routines lantr
, lanhe
, lansy
, lange
depending on matrix #flags
determine effective rank either by SVD method or QR-factorization with pivoting
if overwrite_a is true, a will be overriden in process of calculation
QR method is faster, but could fail to determine rank in some cases
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]]
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]]
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
Solves matrix equation self*x = b
and returns x
Method returns matrix of same size as b
Matrix must be square, number of rows must match b
if overwrite_a is true, a will be overriden in process of calculation
if overwrite_b is true, b will be overriden in process of calculation
Uses LAPACK routines trtrs
, posv
, hesv
, sysv
, gesv
depending on matrix #flags
Calculates singular value decomposition for a matrix
If you call u,s,vt = a.svd
for a matrix m*n a
then
u
will be m*m matrixs
will be Array(T) with size equal to{m, n}.min
vt
will be n*n matrix(u*Mat.diag(a.nrows, a.ncolumns, s)*vt)
will be equal toa
(within calculation tolerance)
if overwrite_a is true, a will be overriden in process of calculation
Uses gesdd
LAPACK routine
Calculates array of singular values for a matrix
if overwrite_a
is true, a will be overriden in process of calculation
Uses gesdd
LAPACK routine
Converts matrix to plain array of elements
if col_major
is true, elements are returned as stored inplace,
otherwise row-major storage is emulated
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]]
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
transposes matrix inplace
Currently, transpose of non-square matrix still allocates temporary buffer
sets element at row i and column j to value, without performing any checks