Module: Numo::Linalg
- Defined in:
- lib/numo/linalg.rb,
lib/numo/linalg/version.rb,
ext/numo/linalg/linalg.c,
ext/numo/linalg/linalg.c
Overview
Numo::Linalg is a subset library from Numo::Linalg consisting only of methods used in Machine Learning algorithms.
Defined Under Namespace
Classes: LapackError
Constant Summary collapse
- VERSION =
The version of numo-linalg-alt you install.
'0.9.1'- OPENBLAS_VERSION =
The version of OpenBLAS used in background library.
rb_str_new_cstr(OPENBLAS_VERSION)
Class Method Summary collapse
-
.blas_char(a, ...) ⇒ String
Returns BLAS char ([sdcz]) defined by data-type of arguments.
-
.cho_fact(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a symmetric / Hermitian positive-definite matrix.
-
.cho_inv(a, uplo: 'U') ⇒ Numo::NArray
Computes the inverse of a symmetric / Hermitian positive-definite matrix using its Cholesky decomposition.
-
.cho_solve(a, b, uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` with the Cholesky factorization of `A`.
-
.cho_solve_banded(ab, b, uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` with the Cholesky factorization of banded matrix `A`.
-
.cholesky(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a symmetric / Hermitian positive-definite matrix.
-
.cholesky_banded(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a banded symmetric / Hermitian positive-definite matrix.
-
.cond(a, ord = nil) ⇒ Numo::NArray
Compute the condition number of a matrix.
-
.coshm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic cosine using the matrix exponential.
-
.cosm(a) ⇒ Numo::NArray
Computes the matrix cosine using the matrix exponential.
-
.det(a) ⇒ Float/Complex
Computes the determinant of matrix.
-
.diagsvd(s, m, n) ⇒ Numo::NArray
Creates a diagonal matrix from the given singular values.
-
.dot(a, b) ⇒ Float|Complex|Numo::NArray
Calculates dot product of two vectors / matrices.
-
.eig(a, left: false, right: true) ⇒ Array<Numo::NArray>
Computes the eigenvalues and right and/or left eigenvectors of a general square matrix.
-
.eig_banded(ab, vals_only: false, vals_range: nil, lower: false) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a symmetric / Hermitian banded matrix.
-
.eigh(a, b = nil, vals_only: false, vals_range: nil, uplo: 'U', turbo: false) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a symmetric / Hermitian matrix by solving an ordinary or generalized eigenvalue problem.
-
.eigh_tridiagonal(d, e, vals_only: false, vals_range: nil) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a real symmetric tridiagonal matrix.
-
.eigvals(a) ⇒ Numo::NArray
Computes the eigenvalues of a general square matrix.
-
.eigvals_banded(ab, vals_range: nil, lower: false) ⇒ Numo::NArray
Computes the eigenvalues a symmetric / Hermitian banded matrix.
-
.eigvalsh(a, b = nil, vals_range: nil, uplo: 'U', turbo: false) ⇒ Numo::NArray
Computes the eigenvalues of a symmetric / Hermitian matrix by solving an ordinary / generalized eigenvalue problem.
-
.eigvalsh_tridiagonal(d, e, vals_range: nil) ⇒ Numo::NArray
Computes the eigenvalues of a real symmetric tridiagonal matrix.
-
.expm(a, ord = 8) ⇒ Numo::NArray
Computes the matrix exponential using a scaling and squaring algorithm with a Pade approximation.
-
.hessenberg(a, calc_q: false) ⇒ Numo::NArray, Array<Numo::NArray, Numo::NArray>
Computes the Hessenberg decomposition of a square matrix.
-
.inv(a, driver: 'getrf', uplo: 'U') ⇒ Numo::NArray
Computes the inverse matrix of a square matrix.
-
.ldl(a, uplo: 'U', hermitian: true) ⇒ Array<Numo::NArray>
Computes the Bunch-Kaufman decomposition of a symmetric / Hermitian matrix.
-
.logm(a) ⇒ Numo::NArray
Computes the matrix logarithm using its eigenvalue decomposition.
-
.lstsq(a, b, driver: 'lsd', rcond: nil) ⇒ Array<Numo::NArray, Float/Complex, Integer, Numo::NArray>
Computes the least-squares solution to a linear matrix equation.
-
.lu(a, permute_l: false) ⇒ Array<Numo::NArray>
Computes the LU decomposition of a matrix using partial pivoting.
-
.lu_fact(a) ⇒ Array<Numo::NArray>
Computes the LU decomposition of a matrix using partial pivoting.
-
.lu_inv(lu, ipiv) ⇒ Numo::NArray
Computes the inverse of a matrix using its LU decomposition.
-
.lu_solve(lu, ipiv, b, trans: 'N') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` using the LU decomposition of `A`.
-
.matmul(a, b) ⇒ Numo::NArray
Computes the matrix multiplication of two arrays.
-
.matrix_balance(a, permute: true, scale: true, separate: false) ⇒ Array<Numo::NArray, Numo::NArray>
Computes a diagonal similarity transformation that balances a square matrix.
-
.matrix_power(a, n) ⇒ Numo::NArray
Computes the matrix ‘a` raised to the power of `n`.
-
.matrix_rank(a, tol: nil, driver: 'svd') ⇒ Integer
Computes the rank of a matrix using SVD.
-
.norm(a, ord = nil, axis: nil, keepdims: false) ⇒ Numo::NArray/Numeric
Computes the matrix or vector norm.
-
.null_space(a, rcond: nil) ⇒ Numo::NArray
Computes an orthonormal basis for the null space of ‘A` using SVD.
-
.orth(a, rcond: nil) ⇒ Numo::NArray
Computes an orthonormal basis for the range of ‘A` using SVD.
-
.orthogonal_procrustes(a, b) ⇒ Array<Numo::NArray, Float>
Computes the orthogonal Procrustes problem.
-
.pinv(a, driver: 'svd', rcond: nil) ⇒ Numo::NArray
Computes the (Moore-Penrose) pseudo-inverse of a matrix using singular value decomposition.
-
.polar(a, side: 'right') ⇒ Array<Numo::NArray>
Computes the polar decomposition of a matrix.
-
.qr(a, mode: 'reduce') ⇒ Numo::NArray+
Computes the QR decomposition of a matrix.
-
.qz(a, b) ⇒ Array<Numo::NArray, Numo::NArray, Numo::NArray, Numo::NArray>
Computes the QZ decomposition (generalized Schur decomposition) of a pair of square matrices.
-
.rq(a, mode: 'full') ⇒ Array<Numo::NArray>/Numo::NArray
Computes the RQ decomposition of a matrix.
-
.schur(a, sort: nil) ⇒ Array<Numo::NArray, Numo::NArray, Integer>
Computes the Schur decomposition of a square matrix.
-
.signm(a) ⇒ Numo::NArray
Computes the matrix sign function using its inverse and square root matrices.
-
.sinhm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic sine using the matrix exponential.
-
.sinm(a) ⇒ Numo::NArray
Computes the matrix sine using the matrix exponential.
-
.slogdet(a) ⇒ Array<Float/Complex>
Computes the sign and natural logarithm of the determinant of a matrix.
-
.solve(a, b, driver: 'gen', uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` from square matrix `A`.
-
.solve_banded(l_u, ab, b) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a banded matrix.
-
.solve_triangular(a, b, lower: false) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a triangular matrix.
-
.solveh_banded(ab, b, lower: false) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a symmetric/Hermitian positive-definite banded matrix.
-
.sqrtm(a) ⇒ Numo::NArray
Computes the square root of a matrix using its eigenvalue decomposition.
-
.svd(a, driver: 'svd', job: 'A') ⇒ Array<Numo::NArray>
Computes the Singular Value Decomposition (SVD) of a matrix: ‘A = U * S * V^T`.
-
.svdvals(a, driver: 'sdd') ⇒ Numo::NArray
Computes the singular values of a matrix.
-
.tanhm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic tangent.
-
.tanm(a) ⇒ Numo::NArray
Computes the matrix tangent.
Class Method Details
.blas_char(a, ...) ⇒ String
Returns BLAS char ([sdcz]) defined by data-type of arguments.
76 77 78 79 80 81 82 83 84 85 86 87 |
# File 'ext/numo/linalg/linalg.c', line 76
static VALUE linalg_blas_char(int argc, VALUE* argv, VALUE self) {
VALUE nary_arr = Qnil;
rb_scan_args(argc, argv, "*", &nary_arr);
const char type = blas_char(nary_arr);
if (type == 'n') {
rb_raise(rb_eTypeError, "invalid data type for BLAS/LAPACK");
return Qnil;
}
return rb_str_new(&type, 1);
}
|
.cho_fact(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a symmetric / Hermitian positive-definite matrix.
1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 |
# File 'lib/numo/linalg.rb', line 1448 def cho_fact(a, uplo: 'U') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] raise ArgumentError, 'uplo must be "U" or "L"' unless %w[U L].include?(uplo) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}potrf" c, info = Numo::Linalg::Lapack.send(fnc, a.dup, uplo: uplo) raise LapackError, "the #{-info}-th argument of #{fnc} had illegal value" if info.negative? if info.positive? raise LapackError, "the leading principal minor of order #{info} is not positive, " \ 'and the factorization could not be completed.' end c end |
.cho_inv(a, uplo: 'U') ⇒ Numo::NArray
Computes the inverse of a symmetric / Hermitian positive-definite matrix using its Cholesky decomposition.
2210 2211 2212 2213 2214 2215 2216 2217 2218 2219 2220 2221 2222 2223 2224 2225 |
# File 'lib/numo/linalg.rb', line 2210 def cho_inv(a, uplo: 'U') bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}potri" inv, info = Numo::Linalg::Lapack.send(fnc, a.dup, uplo: uplo) raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? if info.positive? raise LapackError, "the (#{info - 1}, #{info - 1})-th element of the factor U or L is zero, " \ 'and the inverse could not be computed.' end inv end |
.cho_solve(a, b, uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` with the Cholesky factorization of `A`.
363 364 365 366 367 368 369 370 371 372 373 374 375 376 |
# File 'lib/numo/linalg.rb', line 363 def cho_solve(a, b, uplo: 'U') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] raise Numo::NArray::ShapeError, "incompatible dimensions: a.shape[0] = #{a.shape[0]} != b.shape[0] = #{b.shape[0]}" if a.shape[0] != b.shape[0] bchr = blas_char(a, b) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}potrs" x, info = Numo::Linalg::Lapack.send(fnc, a, b.dup, uplo: uplo) raise LapackError, "the #{-info}-th argument of potrs had illegal value" if info.negative? x end |
.cho_solve_banded(ab, b, uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` with the Cholesky factorization of banded matrix `A`.
399 400 401 402 403 404 405 406 407 408 409 410 411 |
# File 'lib/numo/linalg.rb', line 399 def cho_solve_banded(ab, b, uplo: 'U') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if ab.ndim != 2 raise Numo::NArray::ShapeError, "incompatible dimensions: ab.shape[1] = #{ab.shape[1]} != b.shape[0] = #{b.shape[0]}" if ab.shape[1] != b.shape[0] bchr = blas_char(ab, b) raise ArgumentError, "invalid array type: #{ab.class}" if bchr == 'n' fnc = :"#{bchr}pbtrs" x, info = Numo::Linalg::Lapack.send(fnc, ab, b.dup, uplo: uplo) raise LapackError, "the #{-info}-th argument of potrs had illegal value" if info.negative? x end |
.cholesky(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a symmetric / Hermitian positive-definite matrix.
301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 |
# File 'lib/numo/linalg.rb', line 301 def cholesky(a, uplo: 'U') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}potrf" c, _info = Numo::Linalg::Lapack.send(fnc, a.dup, uplo: uplo) case uplo when 'U' c.triu when 'L' c.tril else raise ArgumentError, "invalid uplo: #{uplo}" end end |
.cholesky_banded(a, uplo: 'U') ⇒ Numo::NArray
Computes the Cholesky decomposition of a banded symmetric / Hermitian positive-definite matrix.
326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 |
# File 'lib/numo/linalg.rb', line 326 def cholesky_banded(a, uplo: 'U') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}pbtrf" c, info = Numo::Linalg::Lapack.send(fnc, a.dup, uplo: uplo) raise LapackError, "the #{-info}-th argument of potrs had illegal value" if info.negative? if info.positive? raise LapackError, "the leading principal minor of order #{info} is not positive, " \ 'and the factorization could not be completed.' end c end |
.cond(a, ord = nil) ⇒ Numo::NArray
Compute the condition number of a matrix.
nil or 2: 2-norm using singular values, ‘fro’: Frobenius norm, ‘info’: infinity norm, and 1: 1-norm.
1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 |
# File 'lib/numo/linalg.rb', line 1906 def cond(a, ord = nil) if ord.nil? || ord == 2 || ord == -2 svals = svdvals(a) if ord == -2 svals[false, -1] / svals[false, 0] else svals[false, 0] / svals[false, -1] end else inv_a = inv(a) norm(a, ord, axis: [-2, -1]) * norm(inv_a, ord, axis: [-2, -1]) end end |
.coshm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic cosine using the matrix exponential.
2125 2126 2127 2128 2129 2130 |
# File 'lib/numo/linalg.rb', line 2125 def coshm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] 0.5 * (expm(a) + expm(-a)) end |
.cosm(a) ⇒ Numo::NArray
Computes the matrix cosine using the matrix exponential.
2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 |
# File 'lib/numo/linalg.rb', line 2079 def cosm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' b = a * 1.0i if %w[z c].include?(bchr) 0.5 * (expm(b) + expm(-b)) else expm(b).real end end |
.det(a) ⇒ Float/Complex
Computes the determinant of matrix.
424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 |
# File 'lib/numo/linalg.rb', line 424 def det(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' getrf = :"#{bchr}getrf" lu, piv, info = Numo::Linalg::Lapack.send(getrf, a.dup) raise LapackError, "the #{-info}-th argument of getrf had illegal value" if info.negative? # info > 0 means the factor U has a zero diagonal element and is singular. # In this case, the determinant is zero. The method should simply return 0.0. # Therefore, the error is not raised here. # raise 'the factor U is singular, ...' if info.positive? det_l = 1 det_u = lu.diagonal.prod det_p = piv.map_with_index { |v, i| v == i + 1 ? 1 : -1 }.prod det_l * det_u * det_p end |
.diagsvd(s, m, n) ⇒ Numo::NArray
Creates a diagonal matrix from the given singular values.
1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 |
# File 'lib/numo/linalg.rb', line 1261 def diagsvd(s, m, n) sz_s = s.size raise ArgumentError, "size of s must be equal to m or n: s.size=#{sz_s}, m=#{m}, n=#{n}" if sz_s != m && sz_s != n mat = s.class.zeros(m, n) if sz_s == m mat[true, 0...m] = s.diag else mat[0...n, true] = s.diag end mat end |
.dot(a, b) ⇒ Float|Complex|Numo::NArray
Calculates dot product of two vectors / matrices.
125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 |
# File 'ext/numo/linalg/linalg.c', line 125
static VALUE linalg_dot(VALUE self, VALUE a_, VALUE b_) {
VALUE a = IsNArray(a_) ? a_ : rb_funcall(numo_cNArray, rb_intern("asarray"), 1, a_);
VALUE b = IsNArray(b_) ? b_ : rb_funcall(numo_cNArray, rb_intern("asarray"), 1, b_);
VALUE arg_arr = rb_ary_new3(2, a, b);
const char type = blas_char(arg_arr);
if (type == 'n') {
rb_raise(rb_eTypeError, "invalid data type for BLAS/LAPACK");
return Qnil;
}
VALUE ret = Qnil;
narray_t* a_nary = NULL;
narray_t* b_nary = NULL;
GetNArray(a, a_nary);
GetNArray(b, b_nary);
const int a_ndim = NA_NDIM(a_nary);
const int b_ndim = NA_NDIM(b_nary);
if (a_ndim == 1) {
if (b_ndim == 1) {
ID fn_id = type == 'c' || type == 'z' ? rb_intern("dotu") : rb_intern("dot");
ret = rb_funcall(rb_mLinalgBlas, rb_intern("call"), 3, ID2SYM(fn_id), a, b);
} else {
VALUE kw_args = rb_hash_new();
if (!RTEST(nary_check_contiguous(b)) &&
RTEST(rb_funcall(b, rb_intern("fortran_contiguous?"), 0))) {
b = rb_funcall(b, rb_intern("transpose"), 0);
rb_hash_aset(kw_args, ID2SYM(rb_intern("trans")), rb_str_new_cstr("N"));
} else {
rb_hash_aset(kw_args, ID2SYM(rb_intern("trans")), rb_str_new_cstr("T"));
}
char fn_name[] = "xgemv";
fn_name[0] = type;
VALUE argv[3] = { b, a, kw_args };
ret = rb_funcallv_kw(rb_mLinalgBlas, rb_intern(fn_name), 3, argv, RB_PASS_KEYWORDS);
}
} else {
if (b_ndim == 1) {
VALUE kw_args = rb_hash_new();
if (!RTEST(nary_check_contiguous(a)) &&
RTEST(rb_funcall(b, rb_intern("fortran_contiguous?"), 0))) {
a = rb_funcall(a, rb_intern("transpose"), 0);
rb_hash_aset(kw_args, ID2SYM(rb_intern("trans")), rb_str_new_cstr("T"));
} else {
rb_hash_aset(kw_args, ID2SYM(rb_intern("trans")), rb_str_new_cstr("N"));
}
char fn_name[] = "xgemv";
fn_name[0] = type;
VALUE argv[3] = { a, b, kw_args };
ret = rb_funcallv_kw(rb_mLinalgBlas, rb_intern(fn_name), 3, argv, RB_PASS_KEYWORDS);
} else {
VALUE kw_args = rb_hash_new();
if (!RTEST(nary_check_contiguous(a)) &&
RTEST(rb_funcall(b, rb_intern("fortran_contiguous?"), 0))) {
a = rb_funcall(a, rb_intern("transpose"), 0);
rb_hash_aset(kw_args, ID2SYM(rb_intern("transa")), rb_str_new_cstr("T"));
} else {
rb_hash_aset(kw_args, ID2SYM(rb_intern("transa")), rb_str_new_cstr("N"));
}
if (!RTEST(nary_check_contiguous(b)) &&
RTEST(rb_funcall(b, rb_intern("fortran_contiguous?"), 0))) {
b = rb_funcall(b, rb_intern("transpose"), 0);
rb_hash_aset(kw_args, ID2SYM(rb_intern("transb")), rb_str_new_cstr("T"));
} else {
rb_hash_aset(kw_args, ID2SYM(rb_intern("transb")), rb_str_new_cstr("N"));
}
char fn_name[] = "xgemm";
fn_name[0] = type;
VALUE argv[3] = { a, b, kw_args };
ret = rb_funcallv_kw(rb_mLinalgBlas, rb_intern(fn_name), 3, argv, RB_PASS_KEYWORDS);
}
}
RB_GC_GUARD(a);
RB_GC_GUARD(b);
return ret;
}
|
.eig(a, left: false, right: true) ⇒ Array<Numo::NArray>
Computes the eigenvalues and right and/or left eigenvectors of a general square matrix.
1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 |
# File 'lib/numo/linalg.rb', line 1620 def eig(a, left: false, right: true) # rubocop:disable Metrics/AbcSize, Metrics/PerceivedComplexity raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise ArgumentError, 'input array a must be square' if a.shape[0] != a.shape[1] jobvl = left ? 'V' : 'N' jobvr = right ? 'V' : 'N' bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}geev" if %w[z c].include?(bchr) w, vl, vr, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvl: jobvl, jobvr: jobvr) else wr, wi, vl, vr, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvl: jobvl, jobvr: jobvr) end raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'the QR algorithm failed to compute all the eigenvalues.' if info.positive? if %w[d s].include?(bchr) w = wr + (wi * 1.0i) ids = wi.gt(0).where unless ids.empty? cast_class = bchr == 'd' ? Numo::DComplex : Numo::SComplex if left tmp = cast_class.cast(vl) tmp[true, ids].imag = vl[true, ids + 1] tmp[true, ids + 1] = tmp[true, ids].conj vl = tmp end if right tmp = cast_class.cast(vr) tmp[true, ids].imag = vr[true, ids + 1] tmp[true, ids + 1] = tmp[true, ids].conj vr = tmp end end end [w, left ? vl : nil, right ? vr : nil] end |
.eig_banded(ab, vals_only: false, vals_range: nil, lower: false) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a symmetric / Hermitian banded matrix.
1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 |
# File 'lib/numo/linalg.rb', line 1738 def eig_banded(ab, vals_only: false, vals_range: nil, lower: false) raise Numo::NArray::ShapeError, 'input array ab must be 2-dimensional' if ab.ndim != 2 bchr = blas_char(ab) raise ArgumentError, "invalid array type: #{ab.class}" if bchr == 'n' jobz = vals_only ? 'N' : 'V' uplo = lower ? 'L' : 'U' fnc = %w[d s].include?(bchr) ? "#{bchr}sbevx" : "#{bchr}hbevx" if vals_range.nil? _, _, vals, vecs, _, info = Numo::Linalg::Lapack.send(fnc.to_sym, ab.dup, range: 'A', jobz: jobz, uplo: uplo) else il = vals_range.first(1)[0] + 1 iu = vals_range.last(1)[0] + 1 _, _, vals, vecs, _, info = Numo::Linalg::Lapack.send(fnc.to_sym, ab.dup, range: 'I', jobz: jobz, uplo: uplo, il: il, iu: iu) end raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? vecs = nil if vals_only [vals, vecs] end |
.eigh(a, b = nil, vals_only: false, vals_range: nil, uplo: 'U', turbo: false) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a symmetric / Hermitian matrix by solving an ordinary or generalized eigenvalue problem.
55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 |
# File 'lib/numo/linalg.rb', line 55 def eigh(a, b = nil, vals_only: false, vals_range: nil, uplo: 'U', turbo: false) # rubocop:disable Metrics/AbcSize, Metrics/CyclomaticComplexity, Metrics/ParameterLists, Metrics/PerceivedComplexity, Lint/UnusedMethodArgument raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] b_given = !b.nil? raise Numo::NArray::ShapeError, 'input array b must be 2-dimensional' if b_given && b.ndim != 2 raise Numo::NArray::ShapeError, 'input array b must be square' if b_given && b.shape[0] != b.shape[1] raise ArgumentError, "invalid array type: #{b.class}" if b_given && blas_char(b) == 'n' bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' jobz = vals_only ? 'N' : 'V' if b_given fnc = %w[d s].include?(bchr) ? "#{bchr}sygv" : "#{bchr}hegv" if vals_range.nil? fnc << 'd' if turbo vecs, _b, vals, _info = Numo::Linalg::Lapack.send(fnc.to_sym, a.dup, b.dup, jobz: jobz) else fnc << 'x' il = vals_range.first(1)[0] + 1 iu = vals_range.last(1)[0] + 1 _a, _b, _m, vals, vecs, _ifail, _info = Numo::Linalg::Lapack.send( fnc.to_sym, a.dup, b.dup, jobz: jobz, range: 'I', il: il, iu: iu ) end else fnc = %w[d s].include?(bchr) ? "#{bchr}syev" : "#{bchr}heev" if vals_range.nil? fnc << 'd' if turbo vecs, vals, _info = Numo::Linalg::Lapack.send(fnc.to_sym, a.dup, jobz: jobz) else fnc << 'r' il = vals_range.first(1)[0] + 1 iu = vals_range.last(1)[0] + 1 _a, _m, vals, vecs, _isuppz, _info = Numo::Linalg::Lapack.send( fnc.to_sym, a.dup, jobz: jobz, range: 'I', il: il, iu: iu ) end end vecs = nil if vals_only [vals, vecs] end |
.eigh_tridiagonal(d, e, vals_only: false, vals_range: nil) ⇒ Array<Numo::NArray>
Computes the eigenvalues and eigenvectors of a real symmetric tridiagonal matrix.
1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 |
# File 'lib/numo/linalg.rb', line 1810 def eigh_tridiagonal(d, e, vals_only: false, vals_range: nil) # rubocop:disable Metrics/AbcSize raise Numo::NArray::ShapeError, 'input array d must be 1-dimensional' if d.ndim != 1 raise Numo::NArray::ShapeError, 'input array e must be 1-dimensional' if e.ndim != 1 if d.shape[0] != e.shape[0] + 1 raise Numo::NArray::ShapeError, "incompatible dimensions: d.shape[0] (#{d.shape[0]}) != e.shape[0] + 1 (#{e.shape[0] + 1})" end if d.is_a?(Numo::DComplex) || d.is_a?(Numo::SComplex) || e.is_a?(Numo::DComplex) || e.is_a?(Numo::SComplex) raise ArgumentError, 'eigh_tridiagonal does not support complex arrays' end bchr = blas_char(d) raise ArgumentError, "invalid array type: #{d.class}" if bchr == 'n' jobz = vals_only ? 'N' : 'V' fnc = "#{bchr}stevx" e_w_zero = e.concatenate(0) if vals_range.nil? _, vals, vecs, _, info = Numo::Linalg::Lapack.send(fnc.to_sym, d.dup, e_w_zero, range: 'A', jobz: jobz) else il = vals_range.first(1)[0] + 1 iu = vals_range.last(1)[0] + 1 _, vals, vecs, _, info = Numo::Linalg::Lapack.send(fnc.to_sym, d.dup, e_w_zero, range: 'I', jobz: jobz, il: il, iu: iu) end raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? vecs = nil if vals_only [vals, vecs] end |
.eigvals(a) ⇒ Numo::NArray
Computes the eigenvalues of a general square matrix.
1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 |
# File 'lib/numo/linalg.rb', line 1667 def eigvals(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}geev" if %w[z c].include?(bchr) w, _vl, _vr, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvl: 'N', jobvr: 'N') else wr, wi, _vl, _vr, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvl: 'N', jobvr: 'N') w = wr + (wi * 1.0i) end raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'the QR algorithm failed to compute all the eigenvalues.' if info.positive? w end |
.eigvals_banded(ab, vals_range: nil, lower: false) ⇒ Numo::NArray
Computes the eigenvalues a symmetric / Hermitian banded matrix.
1772 1773 1774 |
# File 'lib/numo/linalg.rb', line 1772 def eigvals_banded(ab, vals_range: nil, lower: false) eig_banded(ab, vals_only: true, vals_range: vals_range, lower: lower)[0] end |
.eigvalsh(a, b = nil, vals_range: nil, uplo: 'U', turbo: false) ⇒ Numo::NArray
Computes the eigenvalues of a symmetric / Hermitian matrix by solving an ordinary / generalized eigenvalue problem.
1698 1699 1700 |
# File 'lib/numo/linalg.rb', line 1698 def eigvalsh(a, b = nil, vals_range: nil, uplo: 'U', turbo: false) eigh(a, b, vals_only: true, vals_range: vals_range, uplo: uplo, turbo: turbo)[0] end |
.eigvalsh_tridiagonal(d, e, vals_range: nil) ⇒ Numo::NArray
Computes the eigenvalues of a real symmetric tridiagonal matrix.
1855 1856 1857 |
# File 'lib/numo/linalg.rb', line 1855 def eigvalsh_tridiagonal(d, e, vals_range: nil) eigh_tridiagonal(d, e, vals_only: true, vals_range: vals_range)[0] end |
.expm(a, ord = 8) ⇒ Numo::NArray
Computes the matrix exponential using a scaling and squaring algorithm with a Pade approximation.
Reference:
-
Moler and C. Van Loan, “Nineteen Dubious Ways to Compute the Exponential of a Matrix, Twenty-Five Years Later,” SIAM Review, vol. 45, no. 1, pp. 3-49, 2003.
-
2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 |
# File 'lib/numo/linalg.rb', line 2013 def expm(a, ord = 8) # rubocop:disable Metrics/AbcSize raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] norm = a.abs.max n_sqr = norm.positive? ? [0, Math.log2(norm).to_i + 1].max : 0 a /= 2**n_sqr x = a.dup c = 0.5 sgn = 1 nume = a.class.eye(a.shape[0]) + (c * a) deno = a.class.eye(a.shape[0]) - (c * a) (2..ord).each do |k| c *= (ord - k + 1).fdiv(k * ((2 * ord) - k + 1)) x = a.dot(x) c_x = c * x nume += c_x deno += sgn * c_x sgn = -sgn end a_expm = Numo::Linalg.solve(deno, nume) n_sqr.times { a_expm = a_expm.dot(a_expm) } a_expm end |
.hessenberg(a, calc_q: false) ⇒ Numo::NArray, Array<Numo::NArray, Numo::NArray>
Computes the Hessenberg decomposition of a square matrix. The Hessenberg decomposition is given by ‘A = Q * H * Q^H`, where `A` is the input matrix, `Q` is a unitary matrix, and `H` is an upper Hessenberg matrix.
856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 |
# File 'lib/numo/linalg.rb', line 856 def hessenberg(a, calc_q: false) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' func = :"#{bchr}gebal" b, ilo, ihi, _, info = Numo::Linalg::Lapack.send(func, a.dup) raise LapackError, "the #{-info}-th argument of #{func} had illegal value" if info.negative? func = :"#{bchr}gehrd" hq, tau, info = Numo::Linalg::Lapack.send(func, b, ilo: ilo, ihi: ihi) raise LapackError, "the #{-info}-th argument of #{func} had illegal value" if info.negative? h = hq.triu(-1) return h unless calc_q func = %w[d s].include?(bchr) ? :"#{bchr}orghr" : :"#{bchr}unghr" q, info = Numo::Linalg::Lapack.send(func, hq, tau, ilo: ilo, ihi: ihi) raise LapackError, "the #{-info}-th argument of #{func} had illegal value" if info.negative? [h, q] end |
.inv(a, driver: 'getrf', uplo: 'U') ⇒ Numo::NArray
Computes the inverse matrix of a square matrix.
465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 |
# File 'lib/numo/linalg.rb', line 465 def inv(a, driver: 'getrf', uplo: 'U') # rubocop:disable Lint/UnusedMethodArgument raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' getrf = :"#{bchr}getrf" getri = :"#{bchr}getri" lu, piv, info = Numo::Linalg::Lapack.send(getrf, a.dup) raise LapackError, "the #{-info}-th argument of getrf had illegal value" if info.negative? a_inv, info = Numo::Linalg::Lapack.send(getri, lu, piv) raise LapackError, "the #{-info}-th argument of getrf had illegal value" if info.negative? raise LapackError, 'The matrix is singular, and the inverse matrix could not be computed.' if info.positive? a_inv end |
.ldl(a, uplo: 'U', hermitian: true) ⇒ Array<Numo::NArray>
Computes the Bunch-Kaufman decomposition of a symmetric / Hermitian matrix. The factorization has the form ‘A = U * D * U^T` or `A = L * D * L^T`, where `U` (or `L`) is a product of permutation and unit upper (lower) triangular matrices, and `D` is a block diagonal matrix.
1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 |
# File 'lib/numo/linalg.rb', line 1878 def ldl(a, uplo: 'U', hermitian: true) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' complex = bchr =~ /c|z/ fnc = complex && hermitian ? :"#{bchr}hetrf" : :"#{bchr}sytrf" lud = a.dup ipiv, info = Numo::Linalg::Lapack.send(fnc, lud, uplo: uplo) raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? if info.positive? warn("the factorization has been completed, but the D[#{info - 1}, #{info - 1}] is " \ 'exactly zero, indicating that the block diagonal matrix is singular.') end _lud_permutation(lud, ipiv, uplo: uplo, hermitian: hermitian) end |
.logm(a) ⇒ Numo::NArray
Computes the matrix logarithm using its eigenvalue decomposition.
2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 |
# File 'lib/numo/linalg.rb', line 2044 def logm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] ev, vl, = eig(a, left: true, right: false) v = vl.transpose.conj inv_v = Numo::Linalg.inv(v) log_ev = Numo::NMath.log(ev) inv_v.dot(log_ev.diag).dot(v) end |
.lstsq(a, b, driver: 'lsd', rcond: nil) ⇒ Array<Numo::NArray, Float/Complex, Integer, Numo::NArray>
Computes the least-squares solution to a linear matrix equation.
1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 |
# File 'lib/numo/linalg.rb', line 1960 def lstsq(a, b, driver: 'lsd', rcond: nil) # rubocop:disable Lint/UnusedMethodArgument, Metrics/AbcSize raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, "incompatible dimensions: a.shape[0] = #{a.shape[0]} != b.shape[0] = #{b.shape[0]}" if a.shape[0] != b.shape[0] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' m, n = a.shape if m < n if b.ndim == 1 x = Numo::DFloat.zeros(n) x[0...b.size] = b else x = Numo::DFloat.zeros(n, b.shape[1]) x[0...b.shape[0], 0...b.shape[1]] = b end else x = b.dup end fnc = :"#{bchr}gelsd" s, rank, info = Numo::Linalg::Lapack.send(fnc, a.dup, x, rcond: rcond) raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'the algorithm for computing the SVD failed to converge' if info.positive? resids = x.class[] if m > n if rank == n resids = if b.ndim == 1 (x[n..].abs**2).sum(axis: 0) else (x[n..-1, true].abs**2).sum(axis: 0) end end x = if b.ndim == 1 x[false, 0...n] else x[false, 0...n, true] end end [x, resids, rank, s] end |
.lu(a, permute_l: false) ⇒ Array<Numo::NArray>
Computes the LU decomposition of a matrix using partial pivoting.
1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 |
# File 'lib/numo/linalg.rb', line 1367 def lu(a, permute_l: false) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 m, n = a.shape k = [m, n].min lu, piv = lu_fact(a) l = lu.tril.tap { |nary| nary[nary.diag_indices] = 1 }[true, 0...k].dup u = lu.triu[0...k, 0...n].dup perm = a.class.eye(m).tap do |nary| piv.to_a.each_with_index { |i, j| nary[true, [i - 1, j]] = nary[true, [j, i - 1]].dup } end permute_l ? [perm.dot(l), u] : [perm, l, u] end |
.lu_fact(a) ⇒ Array<Numo::NArray>
Computes the LU decomposition of a matrix using partial pivoting.
1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 |
# File 'lib/numo/linalg.rb', line 1386 def lu_fact(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' getrf = :"#{bchr}getrf" lu, piv, info = Numo::Linalg::Lapack.send(getrf, a.dup) raise LapackError, "the #{info.abs}-th argument of getrf had illegal value" if info.negative? if info.positive? warn("the factorization has been completed, but the factor U[#{info - 1}, #{info - 1}] is " \ 'exactly zero, indicating that the matrix is singular.') end [lu, piv] end |
.lu_inv(lu, ipiv) ⇒ Numo::NArray
Computes the inverse of a matrix using its LU decomposition.
2179 2180 2181 2182 2183 2184 2185 2186 2187 2188 2189 2190 |
# File 'lib/numo/linalg.rb', line 2179 def lu_inv(lu, ipiv) bchr = blas_char(lu) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}getri" inv, info = Numo::Linalg::Lapack.send(fnc, lu.dup, ipiv) raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'the matrix is singular and its inverse could not be computed' if info.positive? inv end |
.lu_solve(lu, ipiv, b, trans: 'N') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` using the LU decomposition of `A`.
1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 |
# File 'lib/numo/linalg.rb', line 1426 def lu_solve(lu, ipiv, b, trans: 'N') raise Numo::NArray::ShapeError, 'input array lu must be 2-dimensional' if lu.ndim != 2 raise Numo::NArray::ShapeError, 'input array lu must be square' if lu.shape[0] != lu.shape[1] raise Numo::NArray::ShapeError, "incompatible dimensions: lu.shape[0] = #{lu.shape[0]} != b.shape[0] = #{b.shape[0]}" if lu.shape[0] != b.shape[0] raise ArgumentError, 'trans must be "N", "T", or "C"' unless %w[N T C].include?(trans) bchr = blas_char(lu) raise ArgumentError, "invalid array type: #{lu.class}" if bchr == 'n' getrs = :"#{bchr}getrs" x, info = Numo::Linalg::Lapack.send(getrs, lu, ipiv, b.dup) raise LapackError, "the #{info.abs}-th argument of getrs had illegal value" if info.negative? x end |
.matmul(a, b) ⇒ Numo::NArray
Computes the matrix multiplication of two arrays.
1162 1163 1164 |
# File 'lib/numo/linalg.rb', line 1162 def matmul(a, b) Numo::Linalg::Blas.call(:gemm, a, b) end |
.matrix_balance(a, permute: true, scale: true, separate: false) ⇒ Array<Numo::NArray, Numo::NArray>
Computes a diagonal similarity transformation that balances a square matrix.
1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 |
# File 'lib/numo/linalg.rb', line 1539 def matrix_balance(a, permute: true, scale: true, separate: false) # rubocop:disable Metrics/AbcSize, Metrics/CyclomaticComplexity, Metrics/MethodLength, Metrics/PerceivedComplexity raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 n = a.shape[0] raise ArgumentError, 'input array a must be square' if a.shape[1] != n bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' job = if permute && scale 'B' elsif permute && !scale 'P' elsif !permute && scale 'S' else 'N' end fnc = :"#{bchr}gebal" b, lo, hi, prm_scl, info = Numo::Linalg::Lapack.send(fnc, a.dup, job: job) raise LapackError, "the #{info.abs}-th argument of #{fnc} had illegal value" if info.negative? # convert from Fortran style index to Ruby style index. lo -= 1 hi -= 1 iprm_scl = Numo::Int32.cast(prm_scl) - 1 # extract scaling factors scaler = prm_scl.class.ones(n) scaler[lo...(hi + 1)] = prm_scl[lo...(hi + 1)] # extract permutation indices perm = Numo::Int32.new(n).seq if hi < n - 1 iprm_scl[(hi + 1)...n].to_a.reverse.each.with_index(1) do |s, i| j = n - i next if s == j tmp_ls, tmp_lj = perm[[s, j]].to_a tmp_rj, tmp_rs = perm[[j, s]].to_a perm[[s, j]] = [tmp_rj, tmp_rs] perm[[j, s]] = [tmp_ls, tmp_lj] end end if lo > 0 # rubocop:disable Style/NumericPredicate iprm_scl[0...lo].to_a.each_with_index do |s, j| next if s == j tmp_ls, tmp_lj = perm[[s, j]].to_a tmp_rj, tmp_rs = perm[[j, s]].to_a perm[[s, j]] = [tmp_rj, tmp_rs] perm[[j, s]] = [tmp_ls, tmp_lj] end end return [b, scaler, perm] if separate # construct inverse permutation matrix inv_perm = Numo::Int32.zeros(n) inv_perm[perm] = Numo::Int32.new(n).seq h = scaler.diag[inv_perm, true].dup [b, h] end |
.matrix_power(a, n) ⇒ Numo::NArray
Computes the matrix ‘a` raised to the power of `n`.
1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 |
# File 'lib/numo/linalg.rb', line 1181 def matrix_power(a, n) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] raise ArgumentError, "exponent n must be an integer: #{n}" unless n.is_a?(Integer) if n.zero? a.class.eye(a.shape[0]) elsif n.positive? r = a.dup (n - 1).times { r = Numo::Linalg.matmul(r, a) } r else inv_a = inv(a) r = inv_a.dup (-n - 1).times { r = Numo::Linalg.matmul(r, inv_a) } r end end |
.matrix_rank(a, tol: nil, driver: 'svd') ⇒ Integer
Computes the rank of a matrix using SVD.
1944 1945 1946 1947 1948 1949 1950 |
# File 'lib/numo/linalg.rb', line 1944 def matrix_rank(a, tol: nil, driver: 'svd') return a.ne(0).any? ? 1 : 0 if a.ndim < 2 s = svdvals(a, driver: driver) tol ||= s.max(axis: -1, keepdims: true) * (a.shape[-2..].max * s.class::EPSILON) s.gt(tol).count(axis: -1) end |
.norm(a, ord = nil, axis: nil, keepdims: false) ⇒ Numo::NArray/Numeric
Computes the matrix or vector norm.
| ord | matrix norm | vector norm |
| ----- | ---------------------- | --------------------------- |
| nil | Frobenius norm | 2-norm |
| 'fro' | Frobenius norm | - |
| 'nuc' | nuclear norm | - |
| 'inf' | x.abs.sum(axis:-1).max | x.abs.max |
| 0 | - | (x.ne 0).sum |
| 1 | x.abs.sum(axis:-2).max | same as below |
| 2 | 2-norm (max sing_vals) | same as below |
| other | - | (x.abs**ord).sum**(1.0/ord) |
133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 |
# File 'lib/numo/linalg.rb', line 133 def norm(a, ord = nil, axis: nil, keepdims: false) # rubocop:disable Metrics/AbcSize, Metrics/CyclomaticComplexity, Metrics/MethodLength, Metrics/PerceivedComplexity a = Numo::NArray.asarray(a) unless a.is_a?(Numo::NArray) return 0.0 if a.empty? # for compatibility with Numo::Linalg.norm if ord.is_a?(String) if ord == 'inf' ord = Float::INFINITY elsif ord == '-inf' ord = -Float::INFINITY end end if axis.nil? norm = case a.ndim when 1 Numo::Linalg::Blas.send(:"#{blas_char(a)}nrm2", a) if ord.nil? || ord == 2 when 2 if ord.nil? || ord == 'fro' Numo::Linalg::Lapack.send(:"#{blas_char(a)}lange", a, norm: 'F') elsif ord.is_a?(Numeric) if ord == 1 Numo::Linalg::Lapack.send(:"#{blas_char(a)}lange", a, norm: '1') elsif !ord.infinite?.nil? && ord.infinite?.positive? Numo::Linalg::Lapack.send(:"#{blas_char(a)}lange", a, norm: 'I') end end else if ord.nil? b = a.flatten.dup Numo::Linalg::Blas.send(:"#{blas_char(b)}nrm2", b) end end unless norm.nil? norm = Numo::NArray.asarray(norm).reshape(*([1] * a.ndim)) if keepdims return norm end end if axis.nil? axis = Array.new(a.ndim) { |d| d } else case axis when Integer axis = [axis] when Array, Numo::NArray axis = axis.flatten.to_a else raise ArgumentError, "invalid axis: #{axis}" end end raise ArgumentError, "the number of dimensions of axis is inappropriate for the norm: #{axis.size}" unless [1, 2].include?(axis.size) raise ArgumentError, "axis is out of range: #{axis}" unless axis.all? { |ax| (-a.ndim...a.ndim).cover?(ax) } if axis.size == 1 ord ||= 2 raise ArgumentError, "invalid ord: #{ord}" unless ord.is_a?(Numeric) ord_inf = ord.infinite? if ord_inf.nil? case ord when 0 a.class.cast(a.ne(0)).sum(axis: axis, keepdims: keepdims) when 1 a.abs.sum(axis: axis, keepdims: keepdims) else (a.abs**ord).sum(axis: axis, keepdims: keepdims)**1.fdiv(ord) end elsif ord_inf.positive? a.abs.max(axis: axis, keepdims: keepdims) else a.abs.min(axis: axis, keepdims: keepdims) end else ord ||= 'fro' raise ArgumentError, "invalid ord: #{ord}" unless ord.is_a?(String) || ord.is_a?(Numeric) raise ArgumentError, "invalid axis: #{axis}" if axis.uniq.size == 1 r_axis, c_axis = axis.map { |ax| ax.negative? ? ax + a.ndim : ax } norm = if ord.is_a?(String) raise ArgumentError, "invalid ord: #{ord}" unless %w[fro nuc].include?(ord) if ord == 'fro' Numo::NMath.sqrt((a.abs**2).sum(axis: axis)) else b = a.transpose(c_axis, r_axis).dup gesvd = :"#{blas_char(b)}gesvd" s, = Numo::Linalg::Lapack.send(gesvd, b, jobu: 'N', jobvt: 'N') s.sum(axis: -1) end else ord_inf = ord.infinite? if ord_inf.nil? case ord when -2 b = a.transpose(c_axis, r_axis).dup gesvd = :"#{blas_char(b)}gesvd" s, = Numo::Linalg::Lapack.send(gesvd, b, jobu: 'N', jobvt: 'N') s.min(axis: -1) when -1 c_axis -= 1 if c_axis > r_axis a.abs.sum(axis: r_axis).min(axis: c_axis) when 1 c_axis -= 1 if c_axis > r_axis a.abs.sum(axis: r_axis).max(axis: c_axis) when 2 b = a.transpose(c_axis, r_axis).dup gesvd = :"#{blas_char(b)}gesvd" s, = Numo::Linalg::Lapack.send(gesvd, b, jobu: 'N', jobvt: 'N') s.max(axis: -1) else raise ArgumentError, "invalid ord: #{ord}" end else r_axis -= 1 if r_axis > c_axis if ord_inf.positive? a.abs.sum(axis: c_axis).max(axis: r_axis) else a.abs.sum(axis: c_axis).min(axis: r_axis) end end end if keepdims norm = Numo::NArray.asarray(norm) unless norm.is_a?(Numo::NArray) norm = norm.reshape(*([1] * a.ndim)) end norm end end |
.null_space(a, rcond: nil) ⇒ Numo::NArray
Computes an orthonormal basis for the null space of ‘A` using SVD.
1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 |
# File 'lib/numo/linalg.rb', line 1333 def null_space(a, rcond: nil) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 s, _u, vt = svd(a, driver: 'sdd', job: 'A') tol = if rcond.nil? || rcond.negative? a.shape.max * s.class::EPSILON else rcond end rank = s.gt(tol * s.max).count vt[rank...vt.shape[0], true].conj.transpose.dup end |
.orth(a, rcond: nil) ⇒ Numo::NArray
Computes an orthonormal basis for the range of ‘A` using SVD.
1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 |
# File 'lib/numo/linalg.rb', line 1296 def orth(a, rcond: nil) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 s, u, = svd(a, driver: 'sdd', job: 'S') tol = if rcond.nil? || rcond.negative? a.shape.max * s.class::EPSILON else rcond end rank = s.gt(tol * s.max).count u[true, 0...rank].dup end |
.orthogonal_procrustes(a, b) ⇒ Array<Numo::NArray, Float>
Computes the orthogonal Procrustes problem.
1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 |
# File 'lib/numo/linalg.rb', line 1497 def orthogonal_procrustes(a, b) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array b must be 2-dimensional' if b.ndim != 2 raise Numo::NArray::ShapeError, "incompatible dimensions: a.shape = #{a.shape} != b.shape = #{b.shape}" if a.shape != b.shape m = b.transpose.dot(a.conj).transpose s, u, vt = svd(m, driver: 'svd', job: 'S') r = u.dot(vt) scale = s.sum [r, scale] end |
.pinv(a, driver: 'svd', rcond: nil) ⇒ Numo::NArray
Computes the (Moore-Penrose) pseudo-inverse of a matrix using singular value decomposition.
504 505 506 507 508 509 510 511 |
# File 'lib/numo/linalg.rb', line 504 def pinv(a, driver: 'svd', rcond: nil) s, u, vh = svd(a, driver: driver, job: 'S') rcond = a.shape.max * s.class::EPSILON if rcond.nil? rank = s.gt(rcond * s[0]).count u = u[true, 0...rank] / s[0...rank] u.dot(vh[0...rank, true]).conj.transpose.dup end |
.polar(a, side: 'right') ⇒ Array<Numo::NArray>
Computes the polar decomposition of a matrix.
544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 |
# File 'lib/numo/linalg.rb', line 544 def polar(a, side: 'right') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise ArugumentError, "invalid side: #{side}" unless %w[left right].include?(side) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' s, w, vh = svd(a, driver: 'svd', job: 'S') u = w.dot(vh) p_mat = if side == 'right' vh.transpose.conj.dot(s.diag).dot(vh) else w.dot(s.diag).dot(w.transpose.conj) end [u, p_mat] end |
.qr(a, mode: 'reduce') ⇒ Numo::NArray+
Computes the QR decomposition of a matrix.
597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 |
# File 'lib/numo/linalg.rb', line 597 def qr(a, mode: 'reduce') raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise ArgumentError, "invalid mode: #{mode}" unless %w[reduce r economic raw].include?(mode) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' geqrf = :"#{bchr}geqrf" qr, tau, = Numo::Linalg::Lapack.send(geqrf, a.dup) return [qr, tau] if mode == 'raw' m, n = qr.shape r = m > n && %w[economic raw].include?(mode) ? qr[0...n, true].triu : qr.triu return r if mode == 'r' org_ung_qr = %w[d s].include?(bchr) ? :"#{bchr}orgqr" : :"#{bchr}ungqr" q = if m < n Numo::Linalg::Lapack.send(org_ung_qr, qr[true, 0...m], tau)[0] elsif mode == 'economic' Numo::Linalg::Lapack.send(org_ung_qr, qr, tau)[0] else qqr = a.class.zeros(m, m) qqr[0...m, 0...n] = qr Numo::Linalg::Lapack.send(org_ung_qr, qqr, tau)[0] end [q, r] end |
.qz(a, b) ⇒ Array<Numo::NArray, Numo::NArray, Numo::NArray, Numo::NArray>
Computes the QZ decomposition (generalized Schur decomposition) of a pair of square matrices.
The QZ decomposition is given by ‘A = Q * AA * Z^H` and `B = Q * BB * Z^H`, where `A` and `B` are the input matrices, `Q` and `Z` are unitary matrices, and `AA` and `BB` are upper triangular matrices (or quasi-upper triangular matrices in real case).
727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 |
# File 'lib/numo/linalg.rb', line 727 def qz(a, b) # rubocop:disable Metrics/AbcSize raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array b must be 2-dimensional' if b.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] raise Numo::NArray::ShapeError, 'input array b must be square' if b.shape[0] != b.shape[1] raise Numo::NArray::ShapeError, "incompatible dimensions: a.shape = #{a.shape} != b.shape = #{b.shape}" if a.shape != b.shape bchr = blas_char(a, b) raise ArgumentError, "invalid array type: #{a.class}, #{b.class}" if bchr == 'n' fnc = :"#{bchr}gges" if %w[d s].include?(bchr) aa, bb, _ar, _ai, _beta, q, z, _sdim, info = Numo::Linalg::Lapack.send(fnc, a.dup, b.dup) else aa, bb, _alpha, _beta, q, z, _sdim, info = Numo::Linalg::Lapack.send(fnc, a.dup, b.dup) end n = a.shape[0] raise LapackError, "the #{-info}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'something other than QZ iteration failed.' if info == n + 1 raise LapackError, "reordering failed in #{bchr}tgsen" if info == n + 3 if info == n + 2 raise LapackError, 'after reordering, roundoff changed values of some eigenvalues ' \ 'so that leading eigenvalues in the Generalized Schur form no ' \ 'longer satisfy the sorting condition.' end if info.positive? && info <= n warn('the QZ iteration failed. (a, b) are not in Schur form, ' \ "but alpha[i] and beta[i] for i = #{info},...,n are correct.") end [aa, bb, q, z] end |
.rq(a, mode: 'full') ⇒ Array<Numo::NArray>/Numo::NArray
Computes the RQ decomposition of a matrix.
672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 |
# File 'lib/numo/linalg.rb', line 672 def rq(a, mode: 'full') # rubocop:disable Metrics/AbcSize raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise ArgumentError, "invalid mode: #{mode}" unless %w[full r economic].include?(mode) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}gerqf" rq, tau, info = Numo::Linalg::Lapack.send(fnc, a.dup) raise LapackError, "the #{-info}-th argument of #{fnc} had illegal value" if info.negative? m, n = rq.shape r = rq.triu(n - m).dup r = r[true, (n - m)...n].dup if mode == 'economic' && n > m return r if mode == 'r' fnc = %w[d s].include?(bchr) ? :"#{bchr}orgrq" : :"#{bchr}ungrq" tmp = if n < m rq[(m - n)...m, 0...n].dup elsif mode == 'economic' rq.dup else rq.class.zeros(n, n).tap { |mat| mat[(n - m)...n, true] = rq } end q, info = Numo::Linalg::Lapack.send(fnc, tmp, tau) raise LapackError, "the #{-info}-th argument of #{fnc} had illegal value" if info.negative? [r, q] end |
.schur(a, sort: nil) ⇒ Array<Numo::NArray, Numo::NArray, Integer>
Computes the Schur decomposition of a square matrix. The Schur decomposition is given by ‘A = Z * T * Z^H`, where `A` is the input matrix, `Z` is a unitary matrix, and `T` is an upper triangular matrix (or quasi-upper triangular matrix in real case).
798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 |
# File 'lib/numo/linalg.rb', line 798 def schur(a, sort: nil) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] raise ArgumentError, "invalid sort: #{sort}" unless sort.nil? || %w[lhp rhp iuc ouc].include?(sort) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' fnc = :"#{bchr}gees" if %w[d s].include?(bchr) b, _wr, _wi, v, sdim, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvs: 'V', sort: sort) else b, _w, v, sdim, info = Numo::Linalg::Lapack.send(fnc, a.dup, jobvs: 'V', sort: sort) end n = a.shape[0] raise LapackError, "the #{-info}-th argument of #{fnc} had illegal value" if info.negative? raise LapackError, 'the QR algorithm failed to compute all the eigenvalues.' if info.positive? && info <= n raise LapackError, 'the eigenvalues could not be reordered.' if info == n + 1 if info == n + 2 raise LapackError, 'after reordering, roundoff changed values of some eigenvalues ' \ 'so that leading eigenvalues in the Schur form no longer satisfy ' \ 'the sorting condition.' end [b, v, sdim] end |
.signm(a) ⇒ Numo::NArray
Computes the matrix sign function using its inverse and square root matrices.
2165 2166 2167 2168 2169 2170 2171 2172 |
# File 'lib/numo/linalg.rb', line 2165 def signm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] a_sqrt = sqrtm(a.dot(a)) a_inv = Numo::Linalg.inv(a) a_inv.dot(a_sqrt) end |
.sinhm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic sine using the matrix exponential.
2114 2115 2116 2117 2118 2119 |
# File 'lib/numo/linalg.rb', line 2114 def sinhm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] 0.5 * (expm(a) - expm(-a)) end |
.sinm(a) ⇒ Numo::NArray
Computes the matrix sine using the matrix exponential.
2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 |
# File 'lib/numo/linalg.rb', line 2060 def sinm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' b = a * 1.0i if %w[z c].include?(bchr) -0.5i * (expm(b) - expm(-b)) else expm(b).imag end end |
.slogdet(a) ⇒ Array<Float/Complex>
Computes the sign and natural logarithm of the determinant of a matrix.
1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 |
# File 'lib/numo/linalg.rb', line 1924 def slogdet(a) lu, ipiv = lu_fact(a) dg = lu.diagonal return 0, (-Float::INFINITY) if dg.eq(0).any? idx = ipiv.class.new(ipiv.shape[-1]).seq(1) n_nonzero = ipiv.ne(idx).count(axis: -1) sign = ((-1.0)**(n_nonzero % 2)) * (dg / dg.abs).prod logdet = Numo::NMath.log(dg.abs).sum(axis: -1) [sign, logdet] end |
.solve(a, b, driver: 'gen', uplo: 'U') ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` from square matrix `A`.
909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 |
# File 'lib/numo/linalg.rb', line 909 def solve(a, b, driver: 'gen', uplo: 'U') # rubocop:disable Lint/UnusedMethodArgument raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a, b) raise ArgumentError, "invalid array type: #{a.class}, #{b.class}" if bchr == 'n' gesv = :"#{bchr}gesv" _lu, x, _ipiv, info = Numo::Linalg::Lapack.send(gesv, a.dup, b.dup) raise LapackError, "the #{-info}-th argument of getrf had illegal value" if info.negative? if info.positive? warn('the factorization has been completed, but the factor is singular, ' \ 'so the solution could not be computed.') end x end |
.solve_banded(l_u, ab, b) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a banded matrix.
1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 |
# File 'lib/numo/linalg.rb', line 1010 def solve_banded(l_u, ab, b) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if ab.ndim != 2 unless l_u.is_a?(Array) && l_u.size == 2 && l_u.all?(Integer) raise ArgumentError, 'l_u must be an array of two integers' end kl, ku = l_u raise Numo::NArray::ShapeError, "ab.shape[0] must be equal to l + u + 1: #{ab.shape[0]} != #{kl} + #{ku} + 1" if ab.shape[0] != kl + ku + 1 bchr = blas_char(ab, b) raise ArgumentError, "invalid array type: #{a.class}, #{b.class}" if bchr == 'n' gbsv = :"#{bchr}gbsv" tmp = ab.class.zeros((2 * kl) + ku + 1, ab.shape[1]) tmp[kl..., true] = ab _, x, _, info = Numo::Linalg::Lapack.send(gbsv, tmp, b.dup, kl: kl, ku: ku) raise LapackError, "wrong value is given to the #{-info}-th argument of #{gbsv} used internally" if info.negative? raise LapackError, 'the leading minor of the matrix is not positive definite' if info.positive? x end |
.solve_triangular(a, b, lower: false) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a triangular matrix.
952 953 954 955 956 957 958 959 960 961 962 963 964 965 |
# File 'lib/numo/linalg.rb', line 952 def solve_triangular(a, b, lower: false) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a, b) raise ArgumentError, "invalid array type: #{a.class}, #{b.class}" if bchr == 'n' trtrs = :"#{bchr}trtrs" uplo = lower ? 'L' : 'U' x, info = Numo::Linalg::Lapack.send(trtrs, a, b.dup, uplo: uplo) raise LapackError, "wrong value is given to the #{info}-th argument of #{trtrs} used internally" if info.negative? x end |
.solveh_banded(ab, b, lower: false) ⇒ Numo::NArray
Solves linear equation ‘A * x = b` or `A * X = B` for `x` assuming `A` is a symmetric/Hermitian positive-definite banded matrix.
1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 |
# File 'lib/numo/linalg.rb', line 1073 def solveh_banded(ab, b, lower: false) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if ab.ndim != 2 bchr = blas_char(ab, b) raise ArgumentError, "invalid array type: #{a.class}, #{b.class}" if bchr == 'n' pbsv = :"#{bchr}pbsv" uplo = lower ? 'L' : 'U' x, info = Numo::Linalg::Lapack.send(pbsv, ab.dup, b.dup, uplo: uplo) raise LapackError, "wrong value is given to the #{info}-th argument of #{pbsv} used internally" if info.negative? raise LapackError, 'the leading minor of the matrix is not positive definite' if info.positive? x end |
.sqrtm(a) ⇒ Numo::NArray
Computes the square root of a matrix using its eigenvalue decomposition.
2149 2150 2151 2152 2153 2154 2155 2156 2157 2158 2159 |
# File 'lib/numo/linalg.rb', line 2149 def sqrtm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] ev, vl, = eig(a, left: true, right: false) v = vl.transpose.conj inv_v = Numo::Linalg.inv(v) sqrt_ev = Numo::NMath.sqrt(ev) inv_v.dot(sqrt_ev.diag).dot(v) end |
.svd(a, driver: 'svd', job: 'A') ⇒ Array<Numo::NArray>
Computes the Singular Value Decomposition (SVD) of a matrix: ‘A = U * S * V^T`
1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 |
# File 'lib/numo/linalg.rb', line 1122 def svd(a, driver: 'svd', job: 'A') raise ArgumentError, "invalid job: #{job}" unless /^[ASN]/i.match?(job.to_s) bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' case driver.to_s when 'sdd' gesdd = :"#{bchr}gesdd" s, u, vt, info = Numo::Linalg::Lapack.send(gesdd, a.dup, jobz: job) when 'svd' gesvd = :"#{bchr}gesvd" s, u, vt, info = Numo::Linalg::Lapack.send(gesvd, a.dup, jobu: job, jobvt: job) else raise ArgumentError, "invalid driver: #{driver}" end raise LapackError, "the #{info.abs}-th argument had illegal value" if info.negative? raise LapackError, 'the input array has a NAN entry' if info == -4 raise LapackError, 'the did not converge' if info.positive? [s, u, vt] end |
.svdvals(a, driver: 'sdd') ⇒ Numo::NArray
Computes the singular values of a matrix.
1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 |
# File 'lib/numo/linalg.rb', line 1213 def svdvals(a, driver: 'sdd') bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' case driver.to_s when 'sdd' gesdd = :"#{bchr}gesdd" s, _u, _vt, info = Numo::Linalg::Lapack.send(gesdd, a.dup, jobz: 'N') when 'svd' gesvd = :"#{bchr}gesvd" s, _u, _vt, info = Numo::Linalg::Lapack.send(gesvd, a.dup, jobu: 'N', jobvt: 'N') else raise ArgumentError, "invalid driver: #{driver}" end raise LapackError, "the #{info.abs}-th argument had illegal value" if info.negative? raise LapackError, 'the input array has a NAN entry' if info == -4 raise LapackError, 'the decomposition did not converge' if info.positive? s end |
.tanhm(a) ⇒ Numo::NArray
Computes the matrix hyperbolic tangent.
2136 2137 2138 2139 2140 2141 2142 2143 |
# File 'lib/numo/linalg.rb', line 2136 def tanhm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] a_sinh = sinhm(a) a_cosh = coshm(a) a_sinh.dot(Numo::Linalg.inv(a_cosh)) end |
.tanm(a) ⇒ Numo::NArray
Computes the matrix tangent.
2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 |
# File 'lib/numo/linalg.rb', line 2098 def tanm(a) raise Numo::NArray::ShapeError, 'input array a must be 2-dimensional' if a.ndim != 2 raise Numo::NArray::ShapeError, 'input array a must be square' if a.shape[0] != a.shape[1] bchr = blas_char(a) raise ArgumentError, "invalid array type: #{a.class}" if bchr == 'n' a_sin = sinm(a) a_cos = cosm(a) a_sin.dot(Numo::Linalg.inv(a_cos)) end |