abstract class
LA::Matrix(T)
- LA::Matrix(T)
- Reference
- Object
Overview
class that provide all utility matrix functions
Included Modules
- Enumerable(T)
- LA::LapackHelper
Direct Known Subclasses
Defined in:
linalg/cholesky.crlinalg/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:308
matrix/special_matrix.cr:342
Class Method Summary
-
.arange(start_val : T, end_val : T, delta = 1.0)
Create row from start_val...end_val with step of delta between
-
.block_diag(*args)
Create a block diagonal matrix from provided matrices
-
.circulant(c)
Construct a Circulant matrix
-
.column(values)
Returns single column matrix with elements taken from array
values -
.companion(a)
Create a companion matrix associated with the polynomial whose coefficients are given in
a - .dft(n, scale : DFTScale = DFTScale::None)
-
.diag(nrows, ncolumns, value : Number | Complex)
Returns diagonal matrix of given size with all diagonal elements equal to
value -
.diag(nrows, ncolumns, values)
Returns diagonal matrix of given size with diagonal elements taken from array
values -
.diag(values)
Returns square diagonal matrix with diagonal elements taken from array
values -
.diag(nrows, ncolumns, &)
Returns diagonal matrix of given size with diagonal elements equal to block value
-
.eye(n)
returns identity matrix of size
n -
.fiedler(values)
Returns a symmetric Fiedler matrix
-
.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, "(", ",", "|", ")") -
.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: ")") -
.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) -
.hadamard(n)
Constructs an n-by-n Hadamard matrix, n must be power of 2
- .hankel(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil)
-
.helmert(n, full = false)
Create an Helmert matrix of order n
-
.hilbert(n)
Create a Hilbert matrix of order n.
-
.identity(n)
returns identity matrix of size
n -
.invhilbert(n)
Compute the inverse of the Hilbert matrix of order n.
-
.invpascal(n, kind : PascalKind = PascalKind::Symmetric)
Returns the inverse of the n x n Pascal matrix
-
.kron(a, b)
Returns kroneker product of matrices
-
.leslie(f, s)
Create a Leslie matrix
-
.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") -
.ones(nrows, ncolumns)
Generate matrix of given size with all elements equal to one
-
.pascal(n, kind : PascalKind = PascalKind::Symmetric)
Returns the n x n Pascal matrix.
-
.rand(nrows, ncolumns, rng : Random | Nil = nil)
Generate matrix of given size with elements randomly distributed from range 0.0..1.0
-
.repmat(a : Matrix(T), nrows, ncolumns)
Alias for
#repmat -
.row(values)
Returns single row matrix with elements taken from array
values -
.toeplitz(column : Indexable | Matrix, row : Indexable | Matrix | Nil = nil)
Create a Toeplitz matrix
-
.tri(nrows, ncolumns, k = 0)
Construct (nrows, ncolumns) matrix filled with ones at and below the kth diagonal.
-
.zeros(nrows, ncolumns)
Generate matrix of given size with all elements equal to zero
Macro Summary
-
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 ofMatrixor its descendants -
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 ofMatrixor its descendants
Instance Method Summary
-
#*(k : Number)
Multiplies matrix to scalar
-
#*(k : Complex)
Multiplies matrix to scalar
-
#*(m : Matrix(T))
Matrix product to given m
-
#**(other : Int)
Raises the square matrix to the integer power
other -
#+(k : Number)
Adds scalar value to every element
-
#+(k : Complex)
Adds scalar value to every element
-
#+(m : Matrix(T))
Returns element-wise sum
-
#-(k : Number | Complex)
Substracts scalar value from every element
-
#-(m : Matrix(T))
Returns element-wise substract
-
#-
Multiplies matrix to -1
-
#/(k : Number | Complex)
Divides matrix to scalar
-
#==(other)
Compare with another matrix
-
#[](i : Int32, j : Int32)
Return element of matrix on row i and column j
-
#[](arows : Range(Int32 | Nil, Int32 | Nil), acolumns : Range(Int32 | Nil, Int32 | Nil))
Return submatrix over given ranges.
-
#[](row : Int32, acolumns : Range(Int32 | Nil, Int32 | Nil))
Return submatrix over given ranges.
-
#[](arows : Range(Int32 | Nil, Int32 | Nil), column : Int32)
Return submatrix over given ranges.
-
#[]=(i : Int32, j : Int32, value)
Assign element of matrix on row i and column j
-
#[]=(arows : Range(Int32, Int32), acolumns : Range(Int32, Int32), value)
Assign value to a given submatrix
-
#[]=(row : Int32, acolumns : Range(Int32, Int32), value)
Assign value to a given submatrix
-
#[]=(nrows : Range(Int32, Int32), column : Int32, value)
Assign value to a given submatrix
-
#add(m : Matrix, *, alpha = 1, beta = 1)
Calculate linear combination with matrix
ma.add(m, alpha: alpha, beta: beta)is equal toalpha * a + beta * k, but faster as only one new matrix is allocated -
#add!(k : Number, m : Matrix)
Perform inplace addition with matrix
mmultiplied to scalark -
#almost_eq(other : Matrix(T), eps)
Approximately compare with
othermatrix -
#almost_eq(other : Matrix(T))
Approximately compare with
othermatrix -
#assume!(flag : MatrixFlags, value : Bool = true)
Directly set or reset matrix
flagwithout check SeeMatrixFlagsfor description of flags -
#balance(*, permute = true, scale = true, separate = false)
Balances a square matrix to improve eigenvalue accuracy.
-
#cat(other : Matrix(T), axis : Axis)
Returns concatenation with another matrix by
Axis::Rows(horizontal) orAxis::Columns(vertical) -
#cholesky(*, lower = false, dont_clean = false)
Computes the Cholesky decomposition of a symmetric/Hermitian positive-definite matrix.
-
#chop(eps = self.tolerance)
Converts complex matrix to real one if all imaginary parts are less then
eps, returnsnilotherwise -
#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 -
#columns
Returns Indexable(SubMatrix(T)) that allows iterating over columns
-
#conjt
alias to
#conjtranspose -
#conjt!
alias to
#conjtranspose! -
#conjtranspose
Returns conjurgate transposed matrix
-
#coshm
Matrix hyperbolic cosine
cosh(A) = (exp(A) + exp(-A))/2 -
#cosm
Matrix cosine
cos(A), computed via complex exponential. -
#det
Calculates the determinant of a square matrix.
-
#detect(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance)
Detect if given
aflagsare true or flase for a matrix with toleranceepsUpdate#flagsproperty SeeMatrixFlagsfor description of matrix flags Returnsselffor method chaining -
#detect?(aflags : MatrixFlags = MatrixFlags::All, eps = tolerance)
Detect if given
aflagsare true or flase for a matrix with toleranceepsUpdate#flagsproperty SeeMatrixFlagsfor description of matrix flags Returns True if all given flags are set -
#diag(offset = 0)
Returns
Indexable(T)that allow iterating over k-th diagonal of matrix -
#each(*, all = false, &)
Yields every element of matrix
-
#each_index(*, all = false, &)
Yields every index
-
#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
-
#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
-
#each_with_index(*, all = false, &)
Yields every element of matrix with corresponding row and column
- #eigs(*, left = false)
- #eigs(*, need_left : Bool, need_right : Bool)
-
#eigs(*, b : self, need_left = false, need_right = false)
See
GeneralMatrix#eigss(*, b : GeneralMatrix(T), need_left : Bool, need_right : Bool, overwrite_a = false, overwrite_b = false). - #eigvals
-
#expm
Matrix exponential
exp(A)SeeGeneralMatrix#expmfor details -
#flags : MatrixFlags
Returns flags of matrix (see
MatrixFlags) -
#hcat(other)
Returns horizontal concatenation with another matrix
-
#hessenberg
Computes the Hessenberg form in-place.
-
#hessenberg(*, calc_q = false)
Computes the Hessenberg form of a matrix.
-
#inspect(io)
Converts matrix to string for inspection
-
#inv
Calculate matrix inversion.
-
#kron(b : Matrix(T))
Returns kroneker product with matrix b
-
#lq
See
GeneralMatrix#lq -
#lq_r
See
GeneralMatrix#lq_r - #lu
-
#lu_factor : LUMatrix(T)
Compute pivoted LU decomposition of a matrix and returns it in a compact form, useful for solving linear equation systems.
-
#map(&block : T -> U) forall U
Returns result of appliyng block to every element (without index)
-
#map!(&)
Yields each element (without index) and replace it with returned value
-
#map_with_index(&block : T -> U) forall U
Returns result of appliyng block to every element with index
-
#map_with_index!(&)
Yields each element with index and replace it with returned value
-
#max(axis : Axis)
Calculate maximum over given
axis -
#min(axis : Axis)
Calculate minimum over given
axis -
#ncolumns : Int32
Returns number of columns in matrix
-
#norm(kind : MatrixNorm = MatrixNorm::Frobenius)
Computes the norm of the matrix.
-
#nrows : Int32
Returns number of rows in matrix
-
#pinv
Calculate Moore–Penrose pseudo-inverse of a matrix.
-
#product(axis : Axis)
Perform product over given
axis -
#ql
See
GeneralMatrix#ql -
#ql_r
See
GeneralMatrix#ql_r -
#qr(*, pivoting = false)
See
GeneralMatrix#qr - #qr_r(*, pivoting = false)
- #qr_raw(*, pivoting = false)
-
#qz(b : self)
See
GeneralMatrix#qz -
#reduce(axis : Axis, initial, &)
Perform
#reducefrominitialvalue over givenaxis -
#repmat(arows, acolumns)
Return matrix repeated
arowstimes by vertical andacolumnstimes by horizontal -
#rows
Returns Indexable(SubMatrix(T)) that allows iterating over rows
-
#rq
See
GeneralMatrix#rq -
#rq_r
See
GeneralMatrix#rq_r -
#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") -
#scale!(k : Number | Complex)
Perform inplace multiplication to scalar
k - #schur
-
#shape
Returns shape of matrix in a form of tuple {nrows, ncolumns}
-
#shape_str
Returns shape of matrix as a string
[5x7] -
#sinhm
Matrix hyperbolic sine
sinh(A) = (exp(A) - exp(-A))/2 -
#sinm
Matrix sine
sin(A), computed via complex exponential. -
#size : Tuple(Int32, Int32)
Returns shape of matrix in a form of tuple {nrows, ncolumns}
-
#solve(b : self)
Solves the linear system A * X = B.
-
#solve(b : GeneralMatrix(T), *, overwrite_b = false)
Solves the linear system A * X = B.
-
#square?
Returns True if matrix is square and False otherwise
-
#sum(axis : Axis)
Perform sum over given
axis -
#svd
Computes the singular value decomposition (SVD) of a matrix.
-
#svdvals
Computes only the singular values of a matrix.
-
#t
alias to
#transpose -
#t!
alias to
#transpose! -
#tanhm
Matrix hyperbolic tangent
tanh(A) = sinhm(A) * coshm(A)⁻¹ -
#tanm
Matrix tangent
tan(A) = sinm(A) * cosm(A)⁻¹ -
#to_custom(io, prefix, columns_separator, rows_separator, postfix)
to_custom(io, "[", ",", "],[", "]") Converts a matrix to string with custom format.
-
#to_custom(prefix, columns_separator, rows_separator, postfix)
to_custom(io, "[", ",", "],[", "]") Converts a matrix to string with custom format.
-
#to_general
Creates general matrix with same content.
-
#to_imag
Converts complex matrix to imaginary part
-
#to_matlab(io)
Converts a matrix to matlab format
-
#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]" -
#to_real
Converts complex matrix to real part
-
#to_s(io)
Converts matrix to string, with linefeeds before and after matrix:
-
#to_unsafe
Returns pointer to underlying data
-
#tolerance
Returns estimated tolerance of equality\inequality
-
#trace
Returns sum of diagonal elements (trace of matrix)
-
#transpose
Returns transposed matrix
-
#tril(k = 0)
Same as tril in scipy - returns lower triangular or trapezoidal part of matrix
-
#tril!(k = 0)
Works like a tril in scipy - remove all elements above k-diagonal
-
#triu(k = 0)
Same as triu in scipy - returns upper triangular or trapezoidal part of matrix
-
#triu!(k = 0)
Works like a triu in scipy - remove all elements below k-diagonal
-
#vcat(other)
Returns vertical concatenation with another matrix
Instance methods inherited from module Enumerable(T)
product(initial : Complex)product(initial : Complex, &) product
Class Method Detail
Create row from start_val...end_val with step of delta between
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]]
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
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]
Returns diagonal matrix of given size with all diagonal elements equal to value
Returns diagonal matrix of given size with diagonal elements taken from array values
Returns square diagonal matrix with diagonal elements taken from array values
Returns diagonal matrix of given size with diagonal elements equal to block value
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]
Converts a string of given format to matrix Example:
str = "( 1,2,3 | 4,5,6 | 7,8,9 )"
LA::GMat.from_custom(str, "(", ",", "|", ")")
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: ")")
Converts a string from matlab format to matrix Example:
str = "[1,2,3; 4,5,6; 7,8,9]"
LA::GMat.from_matlab(str)
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]
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]]
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]
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]
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
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]
Returns kroneker product of matrices
Resulting matrix size is {a.nrows*b.nrows, a.ncolumns*b.ncolumns}
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]
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")
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]
Generate matrix of given size with elements randomly distributed from range 0.0..1.0
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]
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]]
Macro Detail
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
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
Instance Method Detail
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
Raises the square matrix to the integer power other
Implementation taken from https://github.com/Exilor/matrix/
Returns element-wise sum
This method raises if another matrix doesn't have same size
Returns element-wise substract
This method raises if another matrix doesn't have same size
Compare with another matrix
Returns True only if all element are exactly equal.
Use #almost_eq If you need approximate equality
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
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]]
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]]
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]]
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
Assign value to a given submatrix
value can be a scalar or a matrix of same size as affected submatrix
Assign value to a given submatrix
value can be a scalar or a matrix of same size as affected submatrix
Assign value to a given submatrix
value can be a scalar or a matrix of same size as affected submatrix
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
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
Approximately compare with other matrix
Returns true if all elements are within eps from corresponding elements of other matrix
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
Directly set or reset matrix flag without check
See MatrixFlags for description of flags
Balances a square matrix to improve eigenvalue accuracy.
Arguments:
- permute (Bool) : If
true, permutes to isolate eigenvalues. Default:true. - scale (Bool) : If
true, scales to improve conditioning. Default:true. - separate (Bool) : If
true, returns scaling factors separately. Default:false.
Returns:
- Tuple(GeneralMatrix(T), GeneralMatrix(T)) : Balanced matrix and scaling factors (or diagonal).
LAPACK Routine:
- Uses
gebal(balance matrix).
Returns concatenation with another matrix by Axis::Rows (horizontal) or Axis::Columns (vertical)
Computes the Cholesky decomposition of a symmetric/Hermitian positive-definite matrix.
Converts the matrix to GeneralMatrix and performs in-place decomposition.
Arguments:
- lower (Bool) : If
true, returns lower triangular factor (L). Otherwise, upper (U). Default:false. - dont_clean (Bool) : If
true, does not zero out the opposite triangle. Default:false.
Returns:
- GeneralMatrix(T) : The lower (L) or upper (U) triangular factor such that A = L*L' or A = U'*U.
Side Effect:
- Marks the original matrix as
PositiveDefiniteif decomposition succeeds.
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
Reset matrix flags to None
Usually is done automatically,
but this method could be needed if internal content was changed using #to_unsafe
Returns conjurgate transposed matrix
result is same as #transpose for real matrices
Calculates the determinant of a square matrix.
Returns:
- T : Determinant value.
See GeneralMatrix#det for details.
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
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
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) }
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 }
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] }
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
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
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 }
See GeneralMatrix#eigs(*, need_left : Bool, need_right : Bool, overwrite_a = false).
This version returns eigenvalues and either left or right eigenvectors (depending on left)
See GeneralMatrix#eigss(*, b : GeneralMatrix(T), need_left : Bool, need_right : Bool, overwrite_a = false, overwrite_b = false).
See GeneralMatrix#eigs(*, need_left : Bool, need_right : Bool, overwrite_a = false).
This version returns only eigenvalues
Computes the Hessenberg form in-place.
Returns:
- GeneralMatrix(T) : Hessenberg matrix.
Computes the Hessenberg form of a matrix.
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.
See GeneralMatrix#hessenberg! for details.
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]
Calculate matrix inversion.
Returns:
- GeneralMatrix(T) : The inverted matrix.
See GeneralMatrix#inv! for details on algorithm.
Compute pivoted LU decomposition of a matrix and returns it in a compact form, useful for solving linear equation systems.
Uses getrf LAPACK routine
Returns result of appliyng block to every element (without index)
Returns result of appliyng block to every element with index
Computes the norm of the matrix.
Arguments:
- kind (MatrixNorm) : Type of norm to compute. Default:
MatrixNorm::Frobenius.
Returns:
- Float : Norm value.
Calculate Moore–Penrose pseudo-inverse of a matrix.
https://en.wikipedia.org/wiki/Moore%E2%80%93Penrose_inverse
Implemented using SVD decomposition.
Returns:
- GeneralMatrix(T) : Pseudo-inverse matrix.
Perform #reduce from initial value over given axis
Return matrix repeated arows times by vertical and acolumns times by horizontal
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")
Solves the linear system A * X = B.
Arguments:
- b (Matrix(T)) : Right-hand side matrix.
Returns:
- GeneralMatrix(T) : Solution matrix X.
Solves the linear system A * X = B.
Arguments:
- b (GeneralMatrix(T)) : Right-hand side matrix.
- overwrite_b (Bool) : If
true, allows overwritingb. Default:false.
Returns:
- GeneralMatrix(T) : Solution matrix X.
See GeneralMatrix#solve for details.
Computes the singular value decomposition (SVD) of a matrix.
Returns:
- Tuple(GeneralMatrix(T), Array(T), GeneralMatrix(T)) : U, singular values, V.t.
See GeneralMatrix#svd for details.
Computes only the singular values of a matrix.
Returns:
- Array(T) : Singular values.
See GeneralMatrix#svdvals for details.
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)"
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)"
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]"
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]
Returns pointer to underlying data
Storage format depends of matrix type This method raises at runtime if matrix doesn't have raw pointer
Returns estimated tolerance of equality\inequality
This value is used by default in #almost_eq compare
Same as tril in scipy - returns lower triangular or trapezoidal part of matrix
Returns a matrix with all elements above k-th diagonal zeroed
Same as triu in scipy - returns upper triangular or trapezoidal part of matrix
Returns a matrix with all elements below k-th diagonal zeroed