abstract class LA::Matrix(T)

Overview

class that provide all utility matrix functions

Included Modules

Direct Known Subclasses

Defined in:

linalg/cholesky.cr
linalg/eig.cr
linalg/lapack_helper.cr
linalg/linalg.cr
linalg/lu.cr
linalg/matfun.cr
linalg/mult.cr
linalg/power.cr
linalg/qr.cr
linalg/rq_lq_ql.cr
linalg/schur.cr
matrix/flag_checks.cr
matrix/formatted_reader.cr
matrix/formatted_writer.cr
matrix/iteration.cr
matrix/matrix.cr
matrix/special_matrix.cr:6
matrix/special_matrix.cr:304
matrix/special_matrix.cr:338

Class Method Summary

Macro Summary

Instance Method Summary

Instance methods inherited from module Enumerable(T)

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

Class Method Detail

def self.arange(start_val : T, end_val : T, delta = 1.0) #

Create row from start_val...end_val with step of delta between


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

Create a block diagonal matrix from provided matrices

Given the inputs A, B and C, the output will have these matrices arranged on the diagonal:

Example:

m = Mat.block_diag(a, b, c)

m will have following structure:

[[a, 0, 0],
 [0, b, 0],
 [0, 0, c]]

[View source]
def self.circulant(c) #

Construct a Circulant matrix

c - first column of matrix

Example:

a = circulant([1, 2, 3])
a.to_aa # => [[1, 3, 2],[2, 1, 3],[3, 2, 1]])

Behaviour copied from scipy


[View source]
def self.column(values) #

Returns single column matrix with elements taken from array values


[View source]
def self.companion(a) #

Create a companion matrix associated with the polynomial whose coefficients are given in a

Behaviour copied from scipy

Example:

Mat.companion([1, -10, 31, -30]) # =>
# LA::GeneralMatrix(Float64) (3x3, None):
# [10.0, -31.0, 30.0]
# [1.0, 0.0, 0.0]
# [0.0, 1.0, 0.0]

[View source]
def self.dft(n, scale : DFTScale = DFTScale::None) #

[View source]
def self.diag(nrows, ncolumns, value : Number | Complex) #

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


[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(values) #

Returns square diagonal matrix 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.eye(n) #

returns identity matrix of size n


[View source]
def self.fiedler(values) #

Returns a symmetric Fiedler matrix

Given an sequence of numbers values, Fiedler matrices have the structure f[i, j] = (values[i] - values[j]).abs, and hence zero diagonals and nonnegative entries. A Fiedler matrix has a dominant positive eigenvalue and other eigenvalues are negative.

Behaviour copied from scipy

Example:

Mat.fiedler([1, 4, 12, 45, 77]) # =>
# LA::GeneralMatrix(Float64) (5x5, Symmetric | Hermitian):
# [0.0, 3.0, 11.0, 44.0, 76.0]
# [3.0, 0.0, 8.0, 41.0, 73.0]
# [11.0, 8.0, 0.0, 33.0, 65.0]
# [44.0, 41.0, 33.0, 0.0, 32.0]
# [76.0, 73.0, 65.0, 32.0, 0.0]

[View source]
def self.from_custom(str : String, prefix, columns_separator, rows_separator, postfix) #

Converts a string of given format to matrix Example:

str = "( 1,2,3 | 4,5,6 | 7,8,9 )"
LA::GMat.from_custom(str, "(", ",", "|", ")")

[View source]
def self.from_custom(io, prefix, columns_separator, rows_separator, postfix) #

Converts a string of given format to matrix Example:

str = "( 1,2,3 | 4,5,6 | 7,8,9 )"
LA::GMat.from_custom(IO::Memory.new(str), prefix: "(", columns_separator: ",", rows_separator: "|", postfix: ")")

[View source]
def self.from_matlab(s) #

Converts a string from matlab format to matrix Example:

str = "[1,2,3; 4,5,6; 7,8,9]"
LA::GMat.from_matlab(str)

[View source]
def self.hadamard(n) #

Constructs an n-by-n Hadamard matrix, n must be power of 2

Behaviour copied from scipy

Example:

Mat.hadamard(4) # =>
# LA::GeneralMatrix(Float64) (4x4,  Symmetric | Hermitian):
# [1.0, 1.0, 1.0, 1.0]
# [1.0, -1.0, 1.0, -1.0]
# [1.0, 1.0, -1.0, -1.0]
# [1.0, -1.0, -1.0, 1.0]

[View source]
def self.hankel(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil) #

Create a Hankel matrix The Hankel matrix has constant anti-diagonals, with .column as its first column and .row as its last row. If `row is nil, then row with elements equal to zero is assumed. Behaviour copied from scipy Examples:

a = Mat.hankel([1, 17, 99])
a.to_aa # => [
# [ 1, 17, 99],
# [17, 99,  0],
# [99,  0,  0]]
a = Mat.hankel([1, 2, 3, 4], [4, 7, 7, 8, 9])
a.to_aa # => [
# [1, 2, 3, 4, 7],
# [2, 3, 4, 7, 7],
# [3, 4, 7, 7, 8],
# [4, 7, 7, 8, 9]]

[View source]
def self.helmert(n, full = false) #

Create an Helmert matrix of order n

If full is true the (n x n) matrix will be returned. Otherwise the submatrix that does not include the first row will be returned

Behaviour copied from scipy

Example:

Mat.helmert(5, full: true) # =>
# LA::GeneralMatrix(Float64) (5x5, Orthogonal):
# [0.4472135954999579, 0.4472135954999579, 0.4472135954999579, 0.4472135954999579, 0.4472135954999579]
# [0.7071067811865476, -0.7071067811865476, 0.0, 0.0, 0.0]
# [0.408248290463863, 0.408248290463863, -0.816496580927726, 0.0, 0.0]
# [0.28867513459481287, 0.28867513459481287, 0.28867513459481287, -0.8660254037844386, 0.0]
# [0.22360679774997896, 0.22360679774997896, 0.22360679774997896, 0.22360679774997896, -0.8944271909999159]

[View source]
def self.hilbert(n) #

Create a Hilbert matrix of order n.

Returns the n by n matrix with entries h[i,j] = 1 / (i + j + 1).

Behaviour copied from scipy

Example:

Mat.hilbert(3) # =>
# LA::GeneralMatrix(Float64) (3x3, Symmetric | Hermitian | PositiveDefinite):
# [1.0, 0.5, 0.3333333333333333]
# [0.5, 0.3333333333333333, 0.25]
# [0.3333333333333333, 0.25, 0.2]

[View source]
def self.identity(n) #

returns identity matrix of size n


[View source]
def self.invhilbert(n) #

Compute the inverse of the Hilbert matrix of order n.

The entries in the inverse of a Hilbert matrix are integers.

Behaviour copied from scipy

Example:

Mat.invhilbert(4) # => LA::GeneralMatrix(Float64) (4x4, Symmetric | Hermitian | PositiveDefinite):
# [16.0, -120.0, 240.0, -140.0]
# [-120.0, 1200.0, -2700.0, 1680.0]
# [240.0, -2700.0, 6480.0, -4200.0]
# [-140.0, 1680.0, -4200.0, 2800.0]

Mat.invhilbert(16)[7, 7] # => 4.247509952853739e+19

[View source]
def self.invpascal(n, kind : PascalKind = PascalKind::Symmetric) #

Returns the inverse of the n x n Pascal matrix

The Pascal matrix is a matrix containing the binomial coefficients as its elements.

kind : see PascalKind

Behaviour copied from scipy

Example:

Mat.invpascal(4) # =>
# LA::GeneralMatrix(Float64) (5x5, Symmetric | Hermitian | PositiveDefinite):
# [5.0, -10.0, 10.0, -5.0, 1.0]
# [-10.0, 30.0, -35.0, 19.0, -4.0]
# [10.0, -35.0, 46.0, -27.0, 6.0]
# [-5.0, 19.0, -27.0, 17.0, -4.0]
# [1.0, -4.0, 6.0, -4.0, 1.0]
Mat.invpascal(5, PascalKind::Lower) # =>
# LA::GeneralMatrix(Float64) (5x5, LowerTriangular):
# [1.0, 0.0, 0.0, 0.0, 0.0]
# [-1.0, 1.0, 0.0, 0.0, 0.0]
# [1.0, -2.0, 1.0, 0.0, 0.0]
# [-1.0, 3.0, -3.0, 1.0, 0.0]
# [1.0, -4.0, 6.0, -4.0, 1.0]

[View source]
def self.kron(a, b) #

Returns kroneker product of matrices

Resulting matrix size is {a.nrows*b.nrows, a.ncolumns*b.ncolumns}


[View source]
def self.leslie(f, s) #

Create a Leslie matrix

Given the length n array of fecundity coefficients f and the length n-1 array of survival coefficients s, return the associated Leslie matrix.

Behaviour copied from scipy

Example:

Mat.leslie([0.1, 2.0, 1.0, 0.1], [0.2, 0.8, 0.7]) # =>
# LA::GeneralMatrix(Float64) (4x4, None):
# [0.1, 2.0, 1.0, 0.1]
# [0.2, 0.0, 0.0, 0.0]
# [0.0, 0.8, 0.0, 0.0]
# [0.0, 0.0, 0.7, 0.0]

[View source]
def self.load_csv(filename) #

Loads a matrix from CSV (comma separated values) file Example:

LA::GMat.rand(30, 30).save_csv("./test.csv")
a = LA::GMat.load_csv("./test.csv")

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

Generate matrix of given size with all elements equal to one


[View source]
def self.pascal(n, kind : PascalKind = PascalKind::Symmetric) #

Returns the n x n Pascal matrix.

The Pascal matrix is a matrix containing the binomial coefficients as its elements.

kind : see PascalKind

Behaviour copied from scipy

Example:

Mat.pascal(4) # =>
# LA::GeneralMatrix(Float64) (4x4, Symmetric | Hermitian | PositiveDefinite):
# [1.0, 1.0, 1.0, 1.0]
# [1.0, 2.0, 3.0, 4.0]
# [1.0, 3.0, 6.0, 10.0]
# [1.0, 4.0, 10.0, 20.0]
Mat.pascal(4, PascalKind::Lower) # =>
# LA::GeneralMatrix(Float64) (4x4, LowerTriangular):
# [1.0, 0.0, 0.0, 0.0]
# [1.0, 1.0, 0.0, 0.0]
# [1.0, 2.0, 1.0, 0.0]
# [1.0, 3.0, 3.0, 1.0]

[View source]
def self.rand(nrows, ncolumns, rng = Random::DEFAULT) #

Generate matrix of given size with elements randomly distributed from range 0.0..1.0


[View source]
def self.repmat(a : Matrix(T), nrows, ncolumns) #

Alias for #repmat


[View source]
def self.row(values) #

Returns single row matrix with elements taken from array values


[View source]
def self.toeplitz(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil) #

Create a Toeplitz matrix

.column - first column of matrix

.row - first row of matrix (if nil, it is assumed that row = column.map(&.conj))

Behaviour copied from scipy

Example:

MatComplex.toeplitz([1, 2, 3], [1, 4, 5, 6]) # =>
# LA::GeneralMatrix(Float64) (3x4, None):
# [1.0, 4.0, 5.0, 6.0]
# [2.0, 1.0, 4.0, 5.0]
# [3.0, 2.0, 1.0, 4.0]
MatComplex.toeplitz([1.0, 2 + 3.i, 4 - 1.i]) # =>
# LA::GeneralMatrix(Complex) (3x3, Hermitian):
# [1.0 + 0.0i, 2.0 - 3.0i, 4.0 + 1.0i]
# [2.0 + 3.0i, 1.0 + 0.0i, 2.0 - 3.0i]
# [4.0 - 1.0i, 2.0 + 3.0i, 1.0 + 0.0i]

[View source]
def self.tri(nrows, ncolumns, k = 0) #

Construct (nrows, ncolumns) matrix filled with ones at and below the kth diagonal.

The matrix has A[i,j] == 1 for j <= i + k

k - Number of subdiagonal below which matrix is filled with ones. k = 0 is the main diagonal, k < 0 subdiagonal and k > 0 superdiagonal.

Behaviour copied from scipy

Example:

Mat.tri(3, 5, 2).to_aa # => [[
# [1, 1, 1, 0, 0],
# [1, 1, 1, 1, 0],
# [1, 1, 1, 1, 1]]
Mat.tri(3, 5, -1).to_aa # => [[
# [0, 0, 0, 0, 0],
# [1, 0, 0, 0, 0],
# [1, 1, 0, 0, 0]]

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

Generate matrix of given size with all elements equal to zero


[View source]

Macro Detail

macro lapack(name, *args, worksize = nil) #

Complex utility macros that simplifies calling certain LAPACK functions It substitute first letter, allocate workareas, raise exception if return value is negative Check source to see supported functions Example:

# Calling *geev to calculate eigenvalues
lapack(geev, 'N'.ord.to_u8, 'N'.ord.to_u8, nrows, a, nrows,
  vals.to_unsafe.as(LibCBLAS::ComplexDouble*),
  Pointer(LibCBLAS::ComplexDouble).null, nrows,
  Pointer(LibCBLAS::ComplexDouble).null, nrows, worksize: [2*nrows])

Note that it should be called only from methods of Matrix or its descendants


[View source]
macro lapack_util(name, worksize, *args) #

Utility macros that simplifies calling certain LAPACK functions Arguments that point to matrix should be passed as matrix(arg) Example:

# Calling *lange to calculate infinity-norm
lapack_util(lange, worksize, 'I', m.nrows, m.ncolumns, matrix(m), m.nrows)

Note that it should be called only from methods of Matrix or its descendants


[View source]

Instance Method Detail

def *(k : Number) #

Multiplies matrix to scalar


[View source]
def *(k : Complex) #

Multiplies matrix to scalar


[View source]
def *(m : Matrix(T)) #

Matrix product to given m

Raises ArgumentError if inner dimensions do not match

This method automatically calls optimal function depending on MatrixFlags.

If one of the matrix is square and triangular - trmm is called

If one of the matrix is symmetric\hermitian - symm/hemm is called

Otherwise - gemm is called


[View source]
def **(other : Int) #

Raises the square matrix to the integer power other

Implementation taken from https://github.com/Exilor/matrix/


[View source]
def +(k : Number) #

Adds scalar value to every element


[View source]
def +(k : Complex) #

Adds scalar value to every element


[View source]
def +(m : Matrix(T)) #

Returns element-wise sum

This method raises if another matrix doesn't have same size


[View source]
def -(k : Number | Complex) #

Substracts scalar value from every element


[View source]
def -(m : Matrix(T)) #

Returns element-wise substract

This method raises if another matrix doesn't have same size


[View source]
def - #

Multiplies matrix to -1


[View source]
def /(k : Number | Complex) #

Divides matrix to scalar


[View source]
def ==(other) #

Compare with another matrix

Returns True only if all element are exactly equal. Use #almost_eq If you need approximate equality


[View source]
def [](i : Int32, j : Int32) #

Return element of matrix on row i and column j

If i or j negative they are counted from end of matrix If i>=nrows or j>=ncolumns exception is raised Use #unsafe_fetch if you need to skip these checks


[View source]
def [](arows : Range(Int32 | Nil, Int32 | Nil), acolumns : Range(Int32 | Nil, Int32 | Nil)) #

Return submatrix over given ranges.

See SubMatrix(T) for details on submatrices Example:

m = GMat32.new([[1, 2, 3, 4], [5, 6, 7, 8], [9, 10, 11, 12], [13, 14, 15, 16]])
m1 = m[1..3, 2..3] # => [[7.0, 8.0], [11.0, 12.0], [15.0, 16.0]]
m2 = m[..3, 2]     # => [[3.0], [7.0], [11.0], [15.0]]

[View source]
def [](row : Int32, acolumns : Range(Int32 | Nil, Int32 | Nil)) #

Return submatrix over given ranges.

See SubMatrix(T) for details on submatrices Example:

m = GMat32.new([[1, 2, 3, 4], [5, 6, 7, 8], [9, 10, 11, 12], [13, 14, 15, 16]])
m1 = m[1..3, 2..3] # => [[7.0, 8.0], [11.0, 12.0], [15.0, 16.0]]
m2 = m[..3, 2]     # => [[3.0], [7.0], [11.0], [15.0]]

[View source]
def [](arows : Range(Int32 | Nil, Int32 | Nil), column : Int32) #

Return submatrix over given ranges.

See SubMatrix(T) for details on submatrices Example:

m = GMat32.new([[1, 2, 3, 4], [5, 6, 7, 8], [9, 10, 11, 12], [13, 14, 15, 16]])
m1 = m[1..3, 2..3] # => [[7.0, 8.0], [11.0, 12.0], [15.0, 16.0]]
m2 = m[..3, 2]     # => [[3.0], [7.0], [11.0], [15.0]]

[View source]
def []=(i : Int32, j : Int32, value) #

Assign element of matrix on row i and column j

If i or j negative they are counted from end of matrix If i>=nrows or j>=ncolumns exception is raised Use #unsafe_set if you need to skip these checks Note this method reset all matrix flags


[View source]
def []=(arows : Range(Int32, Int32), acolumns : Range(Int32, Int32), value) #

Assign value to a given submatrix

value can be a scalar or a matrix of same size as affected submatrix


[View source]
def []=(row : Int32, acolumns : Range(Int32, Int32), value) #

Assign value to a given submatrix

value can be a scalar or a matrix of same size as affected submatrix


[View source]
def []=(nrows : Range(Int32, Int32), column : Int32, value) #

Assign value to a given submatrix

value can be a scalar or a matrix of same size as affected submatrix


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

Calculate linear combination with matrix m a.add(m, alpha: alpha, beta: beta) is equal to alpha * a + beta * k, but faster as only one new matrix is allocated


[View source]
def add!(k : Number, m : Matrix) #

Perform inplace addition with matrix m multiplied to scalar k

a.add!(k, b) is equal to a = a + k * b, but faster as no new matrix is allocated


[View source]
def almost_eq(other : Matrix(T), eps) #

Approximately compare with other matrix

Returns true if all elements are within eps from corresponding elements of other matrix


[View source]
def almost_eq(other : Matrix(T)) #

Approximately compare with other matrix

Returns true if all elements are within eps from corresponding elements of other matrix

Uses #tolerance as an eps by default


[View source]
def assume!(flag : MatrixFlags, value : Bool = true) #

Directly set or reset matrix flag without check See MatrixFlags for description of flags


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

[View source]
def cat(other : Matrix(T), axis : Axis) #

Returns concatenation with another matrix by Axis::Rows (horizontal) or Axis::Columns (vertical)


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

[View source]
def chop(eps = self.tolerance) #

Converts complex matrix to real one if all imaginary parts are less then eps, returns nil otherwise

Returns just matrix if it is already real


[View source]
def clear_flags #

Reset matrix flags to None Usually is done automatically, but this method could be needed if internal content was changed using #to_unsafe


[View source]
def columns #

Returns Indexable(SubMatrix(T)) that allows iterating over columns


[View source]
def conjt #

alias to #conjtranspose


[View source]
def conjt! #

alias to #conjtranspose!


[View source]
def conjtranspose #

Returns conjurgate transposed matrix

result is same as #transpose for real matrices


[View source]
def coshm #

[View source]
def cosm #

optimization idea for noncomplex matrix is from scipy


[View source]
def det #

[View source]
def detect(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance) #

Detect if given aflags are true or flase for a matrix with tolerance eps Update #flags property See MatrixFlags for description of matrix flags Returns self for method chaining


[View source]
def detect?(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance) #

Detect if given aflags are true or flase for a matrix with tolerance eps Update #flags property See MatrixFlags for description of matrix flags Returns True if all given flags are set


[View source]
def diag(offset = 0) #

Returns Indexable(T) that allow iterating over k-th diagonal of matrix

Example:

m = GMat32.new([[-1, 2, 3, 4],
                [5, -6, 7, 8],
                [9, 10, -11, 12]])
m.diag(0).to_a.should eq [-1, -6, -11]
m.diag(1).to_a.should eq [2, 7, 12]
m.diag(2).to_a.should eq [3, 8]
m.diag(3).to_a.should eq [4]
expect_raises(ArgumentError) { m.diag(4) }
m.diag(-1).to_a.should eq [5, 10]
m.diag(-2).to_a.should eq [9]
expect_raises(ArgumentError) { m.diag(-3) }
expect_raises(ArgumentError) { m.diag(-4) }

[View source]
def each(*, all = false, &) #

Yields every element of matrix

all argument controls whether to yield all or non-empty elements for banded\sparse matrices Example: m.each { |v| raise "negative element found" if v < 0 }


[View source]
def each_index(*, all = false, &) #

Yields every index

all argument controls whether to yield all or non-empty elements for banded\sparse matrices Example: m.each_index { |i, j| m[i, j] = -m[i, j] }


[View source]
def each_lower(*, diagonal = true, all = false, &) #

For every element of matrix that is below main diagonal it yields a block with value and with corresponding row and column

This method is useful for symmetric matrices and similar cases include_diagonal argument controls whether to include elements on main diagonal all argument controls whether to yield all or non-empty elements for banded\sparse matrices


[View source]
def each_upper(*, diagonal = true, all = false, &) #

For every element of matrix that is above main diagonal it yields a block with value and with corresponding row and column

This method is useful for symmetric matrices and similar cases include_diagonal argument controls whether to include elements on main diagonal all argument controls whether to yield all or non-empty elements for banded\sparse matrices


[View source]
def each_with_index(*, all = false, &) #

Yields every element of matrix with corresponding row and column

all argument controls whether to yield all or non-empty elements for banded\sparse matrices Example: m.each_with_index { |v, i, j| m2[i, j] = v }


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

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

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

[View source]
def eigvals #

[View source]
abstract def flags : MatrixFlags #

Returns flags of matrix (see MatrixFlags)


[View source]
def hcat(other) #

Returns horizontal concatenation with another matrix


[View source]
def hessenberg #

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

[View source]
def inspect(io) #

Converts matrix to string for inspection

Output looks like:

GeneralMatrix(Float64) (10x10, MatrixFlags::None)
[1, 2, 3, .... 10]
[11, 12, 13, .... 20]
...
[91, 92, 93, .... 100]

[View source]
def inv #

Calculate matrix inversion See #inv! for details on algorithm


[View source]
def kron(b : Matrix(T)) #

Returns kroneker product with matrix b


[View source]
def lq #

[View source]
def lq_r #

[View source]
def lu #

[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.

Uses getrf LAPACK routine


[View source]
def map(&) #

Returns result of appliyng block to every element (without index)


[View source]
def map!(&) #

Yields each element (without index) and replace it with returned value


[View source]
def map_with_index(&) #

Returns result of appliyng block to every element with index


[View source]
def map_with_index!(&) #

Yields each element with index and replace it with returned value


[View source]
def max(axis : Axis) #

Calculate maximum over given axis


[View source]
def min(axis : Axis) #

Calculate minimum over given axis


[View source]
abstract def ncolumns : Int32 #

Returns number of columns in matrix


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

[View source]
abstract def nrows : Int32 #

Returns number of rows in matrix


[View source]
def pinv #

Calculate Moore–Penrose inverse of a matrix

https://en.wikipedia.org/wiki/Moore%E2%80%93Penrose_inverse Implemented using an svd decomposition


[View source]
def product(axis : Axis) #

Perform product over given axis


[View source]
def ql #

[View source]
def ql_r #

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

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

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

[View source]
def qz(b : self) #

[View source]
def reduce(axis : Axis, initial, &) #

Perform #reduce from initial value over given axis


[View source]
def repmat(arows, acolumns) #

Return matrix repeated arows times by vertical and acolumns times by horizontal


[View source]
def rows #

Returns Indexable(SubMatrix(T)) that allows iterating over rows


[View source]
def rq #

[View source]
def rq_r #

[View source]
def save_csv(filename) #

Save a matrix to CSV (comma separated values) file Example:

LA::GMat.rand(30, 30).save_csv("./test.csv")
a = LA::GMat.load_csv("./test.csv")

[View source]
def scale!(k : Number | Complex) #

Perform inplace multiplication to scalar k


[View source]
def schur #

[View source]
def shape #

Returns shape of matrix in a form of tuple {nrows, ncolumns}


[View source]
def shape_str #

Returns shape of matrix as a string [5x7]


[View source]
def sinhm #

[View source]
def sinm #

[View source]
def size : Tuple(Int32, Int32) #

Returns shape of matrix in a form of tuple {nrows, ncolumns}


[View source]
def solve(b : self) #

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

[View source]
def square? #

Returns True if matrix is square and False otherwise


[View source]
def sum(axis : Axis) #

Perform sum over given axis


[View source]
def svd #

[View source]
def svdvals #

[View source]
def t #

alias to #transpose


[View source]
def t! #

alias to #transpose!


[View source]
def tanhm #

[View source]
def tanm #

[View source]
def to_custom(io, prefix, columns_separator, rows_separator, postfix) #

to_custom(io, "[", ",", "],[", "]") Converts a matrix to string with custom format. Example:

a = LA::GMat[[1, 2, 3], [4, 5, 6], [7, 8, 9]]
str = String.build do |io|
  a.to_custom(io, prefix: "(", columns_separator: ",", rows_separator: "|", postfix: ")")
end
str # => "(1,2,3|4,5,6|7,8,9)"

[View source]
def to_custom(prefix, columns_separator, rows_separator, postfix) #

to_custom(io, "[", ",", "],[", "]") Converts a matrix to string with custom format. Example:

a = LA::GMat[[1, 2, 3], [4, 5, 6], [7, 8, 9]]
str = a.to_custom(prefix: "(", columns_separator: ",", rows_separator: "|", postfix: ")")
str # => "(1,2,3|4,5,6|7,8,9)"

[View source]
def to_general #

Creates general matrix with same content. Useful for banded\sparse matrices


[View source]
def to_imag #

Converts complex matrix to imaginary part


[View source]
def to_matlab(io) #

Converts a matrix to matlab format


[View source]
def to_matlab #

Converts a matrix to matlab format string Example:

LA::GMat[[1, 2, 3], [4, 5, 6], [7, 8, 9]].to_matlab # => "[1,2,3; 4,5,6; 7,8,9]"

[View source]
def to_real #

Converts complex matrix to real part


[View source]
def to_s(io) #

Converts matrix to string, with linefeeds before and after matrix:

Output looks like:

[1, 2, 3, .... 10]
[11, 12, 13, .... 20]
...
[91, 92, 93, .... 100]

[View source]
def to_unsafe #

Returns pointer to underlying data

Storage format depends of matrix type This method raises at runtime if matrix doesn't have raw pointer


[View source]
def tolerance #

Returns estimated tolerance of equality\inequality

This value is used by default in #almost_eq compare


[View source]
def trace #

Returns sum of diagonal elements (trace of matrix)


[View source]
def transpose #

Returns transposed matrix


[View source]
def tril(k = 0) #

Same as tril in scipy - returns lower triangular or trapezoidal part of matrix

Returns a matrix with all elements above k-th diagonal zeroed


[View source]
def tril!(k = 0) #

Works like a tril in scipy - remove all elements above k-diagonal


[View source]
def triu(k = 0) #

Same as triu in scipy - returns upper triangular or trapezoidal part of matrix

Returns a matrix with all elements below k-th diagonal zeroed


[View source]
def triu!(k = 0) #

Works like a triu in scipy - remove all elements below k-diagonal


[View source]
def vcat(other) #

Returns vertical concatenation with another matrix


[View source]