Module: Numo::Linalg

Defined in:
ext/numo/linalg/lapack/lapack.c,
ext/numo/linalg/lapack/lapack_c.c,
ext/numo/linalg/lapack/lapack_d.c,
ext/numo/linalg/lapack/lapack_s.c,
ext/numo/linalg/lapack/lapack_z.c,
ext/numo/linalg/blas/blas.c,
ext/numo/linalg/blas/blas_c.c,
ext/numo/linalg/blas/blas_d.c,
ext/numo/linalg/blas/blas_s.c,
ext/numo/linalg/blas/blas_z.c,
lib/numo/linalg/function.rb,
lib/numo/linalg/loader.rb,
lib/numo/linalg/version.rb

Defined Under Namespace

Modules: Blas, Lapack, Loader Classes: LapackError

Constant Summary

BLAS_CHAR =
{
 SFloat => "s",
 DFloat => "d",
 SComplex => "c",
 DComplex => "z",
}
VERSION =
"0.1.0"

Class Method Summary collapse

Class Method Details

.blas_char(*args) ⇒ Object



57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
# File 'lib/numo/linalg/function.rb', line 57

def blas_char(*args)
  t = Float
  args.each do |a|
    k =
      case a
      when NArray
        a.class
      when Array
        NArray.array_type(a)
      end
    if k && k < NArray
      t = k::UPCAST[t]
    end
  end
  BLAS_CHAR[t] || raise(TypeError,"invalid data type for BLAS/LAPACK")
end

.cho_fact(a, uplo: 'U') ⇒ Numo::NArray

Computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix A. The factorization has the form

A = U**H * U,  if UPLO = 'U', or
A = L  * L**H,  if UPLO = 'L',

where U is an upper triangular matrix and L is lower triangular

Parameters:

  • a (Numo::NArray)

    n-by-n symmetric matrix A (>= 2-dimensinal NArray)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle.

Returns:

  • (Numo::NArray)

    the factor U or L.



418
419
420
# File 'lib/numo/linalg/function.rb', line 418

def cho_fact(a, uplo:'U')
  Lapack.call(:potrf, a, uplo:uplo)[0]
end

.cho_inv(a, uplo: 'U') ⇒ Numo::NArray

Computes the inverse of a symmetric/Hermitian positive definite matrix A using the Cholesky factorization A = U**T*U or A = L*L**T computed by Linalg.cho_fact.

Parameters:

  • a (Numo::NArray)

    the triangular factor U or L from the Cholesky factorization, as computed by Linalg.cho_fact.

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle.

Returns:

  • (Numo::NArray)

    the upper or lower triangle of the (symmetric) inverse of A.



434
435
436
# File 'lib/numo/linalg/function.rb', line 434

def cho_inv(a, uplo:'U')
  Lapack.call(:potri, a, uplo:uplo)[0]
end

.cho_solve(a, b, uplo: 'U') ⇒ Numo::NArray

Solves a system of linear equations A*X = B with a symmetric/Hermitian positive definite matrix A using the Cholesky factorization A = U**T*U or A = L*L**T computed by Linalg.cho_fact.

Parameters:

  • a (Numo::NArray)

    the triangular factor U or L from the Cholesky factorization, as computed by Linalg.cho_fact.

  • b (Numo::NArray)

    the right hand side matrix B.

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle.

Returns:

  • (Numo::NArray)

    the solution matrix X.



450
451
452
# File 'lib/numo/linalg/function.rb', line 450

def cho_solve(a, b, uplo:'U')
  Lapack.call(:potrs, a, b, uplo:uplo)[0]
end

.cond(a, ord = nil) ⇒ Numo::NArray

Compute the condition number of a matrix using the norm with one of the following order.

|  ord  |  matrix norm           |
| ----- | ---------------------- |
|  nil  | 2-norm using SVD       |
| 'fro' | Frobenius norm         |
| 'inf' | x.abs.sum(axis:-1).max |
|    1  | x.abs.sum(axis:-2).max |
|    2  | 2-norm (max sing_vals) |

Examples:

a = Numo::DFloat[[1, 0, -1], [0, 1, 0], [1, 0, 1]]
=> Numo::DFloat#shape=[3,3]
[[1, 0, -1],
 [0, 1, 0],
 [1, 0, 1]]
LA = Numo::Linalg
LA.cond(a)
=> 1.4142135623730951
LA.cond(a, 'fro')
=> 3.1622776601683795
LA.cond(a, 'inf')
=> 2.0
LA.cond(a, '-inf')
=> 1.0
LA.cond(a, 1)
=> 2.0
LA.cond(a, -1)
=> 1.0
LA.cond(a, 2)
=> 1.4142135623730951
LA.cond(a, -2)
=> 0.7071067811865475
(LA.svdvals(a)).min*(LA.svdvals(LA.inv(a))).min
=> 0.7071067811865475

Parameters:

  • a (Numo::NArray)

    matrix or vector (>= 1-dimensinal NArray)

  • ord (String or Symbol) (defaults to: nil)

    Order of the norm.

Returns:

  • (Numo::NArray)

    cond result



711
712
713
714
715
716
717
718
# File 'lib/numo/linalg/function.rb', line 711

def cond(a,ord=nil)
  if ord.nil?
    s = svdvals(a)
    s[false, 0]/s[false, -1]
  else
    norm(a, ord, axis:[-2,-1]) * norm(inv(a), ord, axis:[-2,-1])
  end
end

.det(a) ⇒ Float or Complex or Numo::NArray

Determinant of a matrix

Parameters:

  • a (Numo::NArray)

    matrix (>= 2-dimensional NArray)

Returns:

  • (Float or Complex or Numo::NArray)


725
726
727
728
729
730
731
# File 'lib/numo/linalg/function.rb', line 725

def det(a)
  lu, piv, = Lapack.call(:getrf, a)
  idx = piv.new_narray.store(piv.class.new(piv.shape[-1]).seq(1))
  m = piv.eq(idx).count_false(axis:-1) % 2
  sign = m * -2 + 1
  lu.diagonal.prod(axis:-1) * sign
end

.dot(a, b) ⇒ Numo::NArray

Dot product.

Parameters:

  • a (Numo::NArray)

    matrix or vector (>= 1-dimensinal NArray)

  • b (Numo::NArray)

    matrix or vector (>= 1-dimensinal NArray)

Returns:

  • (Numo::NArray)

    result of dot product



82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
# File 'lib/numo/linalg/function.rb', line 82

def dot(a, b)
  a = NArray.asarray(a)
  b = NArray.asarray(b)
  case a.ndim
  when 1
    case b.ndim
    when 1
      Blas.call(:dot, a, b)
    else
      Blas.call(:gemv, b, a, trans:'t')
    end
  else
    case b.ndim
    when 1
      Blas.call(:gemv, a, b)
    else
      Blas.call(:gemm, a, b)
    end
  end
end

.eig(a, left: false, right: true) ⇒ [w,vl,vr]

Computes the eigenvalues and, optionally, the left and/or right eigenvectors for a square nonsymmetric matrix A.

Parameters:

  • a (Numo::NArray)

    square nonsymmetric matrix (>= 2-dimensinal NArray)

  • left (Bool)

    (optional) If true, left eigenvectors are computed.

  • right (Bool)

    (optional) If true, right eigenvectors are computed.

Returns:

  • ([w,vl,vr])
    • w [Numo::NArray] – The eigenvalues.
    • vl [Numo::NArray] – The left eigenvectors if left is true, otherwise nil.
    • vr [Numo::NArray] – The right eigenvectors if right is true, otherwise nil.


468
469
470
471
472
473
474
475
476
477
478
479
480
# File 'lib/numo/linalg/function.rb', line 468

def eig(a, left:false, right:true)
  jobvl, jobvr = left, right
  case blas_char(a)
  when /c|z/
    w, vl, vr, info = Lapack.call(:geev, a, jobvl:jobvl, jobvr:jobvr)
  else
    wr, wi, vl, vr, info = Lapack.call(:geev, a, jobvl:jobvl, jobvr:jobvr)
    w  = wr + wi * Complex::I
    vl = _make_complex_eigvecs(w,vl) if left
    vr = _make_complex_eigvecs(w,vr) if right
  end
  [w,vl,vr] #.compact
end

.eigh(a, vals_only: false, uplo: false, turbo: false) ⇒ [w,v]

Computes the eigenvalues and, optionally, the left and/or right eigenvectors for a square symmetric/hermitian matrix A.

Parameters:

  • a (Numo::NArray)

    square nonsymmetric matrix (>= 2-dimensinal NArray)

  • values_only (Bool)

    (optional) If false, eigenvectors are computed.

  • uplo (String or Symbol)

    (optional, default=’U’) Access upper (‘U’) or lower (‘L’) triangle.

Returns:

  • ([w,v])
    • w [Numo::NArray] – The eigenvalues.
    • v [Numo::NArray] – The eigenvectors if vals_only is false, otherwise nil.


493
494
495
496
497
498
499
500
501
502
503
# File 'lib/numo/linalg/function.rb', line 493

def eigh(a, vals_only:false, uplo:false, turbo:false)
  jobz = vals_only ? 'N' : 'V' # jobz: Compute eigenvalues and eigenvectors.
  case blas_char(a)
  when /c|z/
    func = turbo ? :hegv : :heev
  else
    func = turbo ? :sygv : :syev
  end
  w, v, = Lapack.call(func, a, uplo:uplo, jobz:jobz)
  [w,v] #.compact
end

.eigvals(a) ⇒ Numo::NArray

Computes the eigenvalues only for a square nonsymmetric matrix A.

Parameters:

  • a (Numo::NArray)

    square nonsymmetric matrix (>= 2-dimensinal NArray)

Returns:

  • (Numo::NArray)

    eigenvalues



510
511
512
513
514
515
516
517
518
519
520
# File 'lib/numo/linalg/function.rb', line 510

def eigvals(a)
  jobvl, jobvr = 'N','N'
  case blas_char(a)
  when /c|z/
    w, = Lapack.call(:geev, a, jobvl:jobvl, jobvr:jobvr)
  else
    wr, wi, = Lapack.call(:geev, a, jobvl:jobvl, jobvr:jobvr)
    w  = wr + wi * Complex::I
  end
  w
end

.eigvalsh(a, uplo: false, turbo: false) ⇒ Numo::NArray

Computes the eigenvalues for a square symmetric/hermitian matrix A.

Parameters:

  • a (Numo::NArray)

    square symmetric/hermitian matrix (>= 2-dimensinal NArray)

  • uplo (String or Symbol)

    (optional, default=’U’) Access upper (‘U’) or lower (‘L’) triangle.

Returns:

  • (Numo::NArray)

    eigenvalues



530
531
532
533
534
535
536
537
538
539
# File 'lib/numo/linalg/function.rb', line 530

def eigvalsh(a, uplo:false, turbo:false)
  jobz = 'N' # jobz: Compute eigenvalues and eigenvectors.
  case blas_char(a)
  when /c|z/
    func = turbo ? :hegv : :heev
  else
    func = turbo ? :sygv : :syev
  end
  Lapack.call(func, a, uplo:uplo, jobz:jobz)[0]
end

.inv(a, driver: "getrf", uplo: 'U') ⇒ Numo::NArray

Inverse matrix from square matrix a

Examples:

Numo::Linalg.inv(a,driver:'getrf')
=> Numo::DFloat#shape=[2,2]
[[-2, 1],
 [1.5, -0.5]]
a.dot(Numo::Linalg.inv(a,driver:'getrf'))
=> Numo::DFloat#shape=[2,2]
[[1, 0],
 [8.88178e-16, 1]]

Parameters:

  • a (Numo::NArray)

    n-by-n square matrix (>= 2-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK diriver (‘ge’|’sy’|’he’|’po’) + (“sv”|”trf”) (optional, default=’getrf’)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle. (omitted when driver:”ge”)

Returns:

  • (Numo::NArray)

    The inverse matrix.



829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
# File 'lib/numo/linalg/function.rb', line 829

def inv(a, driver:"getrf", uplo:'U')
  case driver
  when /(ge|sy|he|po)sv$/
    d = $1
    b = a.new_zeros.eye
    solve(a, b, driver:d, uplo:uplo)
  when /(ge|sy|he)tr[fi]$/
    d = $1
    lu, piv = lu_fact(a, driver:d, uplo:uplo)
    lu_inv(lu, piv, driver:d, uplo:uplo)
  when /potr[fi]$/
    lu = cho_fact(a, uplo:uplo)
    cho_inv(lu, uplo:uplo)
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.lstsq(a, b, driver: 'lsd', rcond: -1)) ⇒ [x, resids, rank, s]

Computes the minimum-norm solution to a linear least squares problem:

    minimize 2-norm(| b - A*x |)

using the singular value decomposition (SVD) of A. A is an M-by-N matrix which may be rank-deficient.

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • b (Numo::NArray)

    m-by-nrhs right-hand-side matrix b (>= 1-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK driver from ‘lsd’,’lss’,’lsy’ (optional, default=’lsd’)

  • rcond (Float)

    (optional, default=-1) RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

Returns:

  • ([x, resids, rank, s])
    • x – The solution matrix/vector X.
    • resids – Sums of residues, squared 2-norm for each column in b - a x. If matrix_rank(a) < N or > M, or ‘gelsy’ is used, this is an empty array.
    • rank – The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).
    • s – The singular values of A in decreasing order. Returns nil if ‘gelsy’ is used.


873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
# File 'lib/numo/linalg/function.rb', line 873

def lstsq(a, b, driver:'lsd', rcond:-1)
  a = NArray.asarray(a)
  b = NArray.asarray(b)
  b_orig = nil
  if b.shape.size==1
    b_orig = b
    b = b_orig[true,:new]
  end
  m = a.shape[-2]
  n = a.shape[-1]
  #nrhs = b.shape[-1]
  if m != b.shape[-2]
    raise NArray::ShapeError, "size mismatch: A-row and B-row"
  end
  if m < n   # need to extend b matrix
    shp = b.shape
    shp[-2] = n
    b2 = b.class.zeros(*shp)
    b2[false,0...m,true] = b
    b = b2
  end
  case driver.to_s
  when /^(ge)?lsd$/i
    # x, s, rank, info
    x, s, rank, = Lapack.call(:gelsd, a, b, rcond:rcond)
  when /^(ge)?lss$/i
    # v, x, s, rank, info
    _, x, s, rank, = Lapack.call(:gelss, a, b, rcond:rcond)
  when /^(ge)?lsy$/i
    jpvt = Int32.zeros(*a[false,0,true].shape)
    # v, x, jpvt, rank, info
    _, x, _, rank, = Lapack.call(:gelsy, a, b, jpvt, rcond:rcond)
    s = nil
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
  resids = nil
  if m > n
    if /ls(d|s)$/i =~ driver
      case rank
      when n
        resids = (x[n..-1,true].abs**2).sum(axis:0)
      when NArray
        if true
          resids = (x[false,n..-1,true].abs**2).sum(axis:-2)
        else
          resids = x[false,0,true].new_zeros
          mask = rank.eq(n)
          # NArray does not suppurt this yet.
          resids[mask,true] = (x[mask,n..-1,true].abs**2).sum(axis:-2)
        end
      end
    end
    x = x[false,0...n,true]
  end
  if b_orig && b_orig.shape.size==1
    x = x[true,0]
    resids &&= resids[false,0]
  end
  [x, resids, rank, s]
end

.lu_fact(a, driver: "gen", uplo: "U") ⇒ [lu, ipiv]

Computes an LU factorization of a M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form

A = P * L * U

where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK diriver from ‘gen’,’sym’,’her’. (optional, default=’gen’)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle. (omitted when driver:”gen”)

Returns:

  • ([lu, ipiv])
    • lu [Numo::NArray] – The factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
    • ipiv [Numo::NArray] – The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).


321
322
323
324
325
326
327
328
329
330
331
# File 'lib/numo/linalg/function.rb', line 321

def lu_fact(a, driver:"gen", uplo:"U")
  case driver.to_s
  when /^gen?(trf)?$/i
    Lapack.call(:getrf, a)[0..1]
  when /^(sym?|her?)(trf)?$/i
    func = driver[0..2].downcase+"trf"
    Lapack.call(func, a, uplo:uplo)[0..1]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.lu_inv(lu, ipiv, driver: "gen", uplo: "U") ⇒ Numo::NArray

Computes the inverse of a matrix using the LU factorization computed by Numo::Linalg.lu_fact.

This method inverts U and then computes inv(A) by solving the system

inv(A)*L = inv(U)

for inv(A).

Parameters:

  • lu (Numo::NArray)

    matrix containing the factors L and U from the factorization A = P*L*U as computed by Numo::Linalg.lu_fact.

  • ipiv (Numo::NArray)

    The pivot indices from Numo::Linalg.lu_fact; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).

  • driver (String or Symbol)

    choose LAPACK diriver from ‘gen’,’sym’,’her’. (optional, default=’gen’)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle. (omitted when driver:”gen”)

Returns:

  • (Numo::NArray)

    the inverse of the original matrix A.



354
355
356
357
358
359
360
361
362
363
364
# File 'lib/numo/linalg/function.rb', line 354

def lu_inv(lu, ipiv, driver:"gen", uplo:"U")
  case driver.to_s
  when /^gen?(tri)?$/i
    Lapack.call(:getri, lu, ipiv)[0]
  when /^(sym?|her?)(tri)?$/i
    func = driver[0..2].downcase+"tri"
    Lapack.call(func, lu, ipiv, uplo:uplo)[0]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.lu_solve(lu, ipiv, b, driver: "gen", uplo: "U", trans: "N") ⇒ Numo::NArray

Solves a system of linear equations

A * X = B  or  A**T * X = B

with a N-by-N matrix A using the LU factorization computed by Numo::Linalg.lu_fact

Parameters:

  • lu (Numo::NArray)

    matrix containing the factors L and U from the factorization A = P*L*U as computed by Numo::Linalg.lu_fact.

  • ipiv (Numo::NArray)

    The pivot indices from Numo::Linalg.lu_fact; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).

  • b (Numo::NArray)

    the right hand side matrix B.

  • driver (String or Symbol)

    choose LAPACK diriver from ‘gen’,’sym’,’her’. (optional, default=’gen’)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle. (omitted when driver:”gen”)

  • trans (String or Symbol)

    Specifies the form of the system of equations (omitted if not driver:”gen”):

    • If ‘N’: A * X = B (No transpose).
    • If ‘T’: A*\*T* X = B (Transpose).
    • If ‘C’: A*\*T* X = B (Conjugate transpose = Transpose).

Returns:

  • (Numo::NArray)

    the solution matrix X.



393
394
395
396
397
398
399
400
401
402
403
# File 'lib/numo/linalg/function.rb', line 393

def lu_solve(lu, ipiv, b, driver:"gen", uplo:"U", trans:"N")
  case driver.to_s
  when /^gen?(trs)?$/i
    Lapack.call(:getrs, lu, ipiv, b, trans:trans)[0]
  when /^(sym?|her?)(trs)?$/i
    func = driver[0..2].downcase+"trs"
    Lapack.call(func, lu, ipiv, b, uplo:uplo)[0]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.matmul(a, b) ⇒ Numo::NArray

Matrix product.

Parameters:

  • a (Numo::NArray)

    matrix (>= 2-dimensinal NArray)

  • b (Numo::NArray)

    matrix (>= 2-dimensinal NArray)

Returns:

  • (Numo::NArray)

    result of matrix product



107
108
109
# File 'lib/numo/linalg/function.rb', line 107

def matmul(a, b)
  Blas.call(:gemm, a, b)
end

.matrix_power(a, n) ⇒ Object

Compute a square matrix a to the power n.

  • If n > 0: return a**n.
  • If n == 0: return identity matrix.
  • If n < 0: return (a*\*-1)*\*n.abs.

Examples:

i = Numo::DFloat[[0, 1], [-1, 0]]
=> Numo::DFloat#shape=[2,2]
[[0, 1],
 [-1, 0]]
Numo::Linalg.matrix_power(i,3)
=> Numo::DFloat#shape=[2,2]
[[0, -1],
 [1, 0]]
Numo::Linalg.matrix_power(i,0)
=> Numo::DFloat#shape=[2,2]
[[1, 0],
 [0, 1]]
Numo::Linalg.matrix_power(i,-3)
=> Numo::DFloat#shape=[2,2]
[[0, 1],
 [-1, 0]]

q = Numo::DFloat.zeros(4,4)
q[0..1,0..1] = -i
q[2..3,2..3] = i
q
=> Numo::DFloat#shape=[4,4]
[[-0, -1, 0, 0],
 [1, -0, 0, 0],
 [0, 0, 0, 1],
 [0, 0, -1, 0]]
Numo::Linalg.matrix_power(q,2)
=> Numo::DFloat#shape=[4,4]
[[-1, 0, 0, 0],
 [0, -1, 0, 0],
 [0, 0, -1, 0],
 [0, 0, 0, -1]]

Parameters:

  • a (Numo::NArray)

    square matrix (>= 2-dimensinal NArray).

  • n (Integer)

    the exponent.



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
# File 'lib/numo/linalg/function.rb', line 153

def matrix_power(a, n)
  a = NArray.asarray(a)
  m,k = a.shape[-2..-1]
  unless m==k
    raise NArray::ShapeError, "input must be a square array"
  end
  unless Integer===n
    raise ArgumentError, "exponent must be an integer"
  end
  if n == 0
    return a.class.eye(m)
  elsif n < 0
    a = inv(a)
    n = n.abs
  end
  if n <= 3
    r = a
    (n-1).times do
      r = matmul(r,a)
    end
  else
    while (n & 1) == 0
      a = matmul(a,a)
      n >>= 1
    end
    r = a
    while n != 0
      a = matmul(a,a)
      n >>= 1
      if (n & 1) != 0
        r = matmul(r,a)
      end
    end
  end
  r
end

.matrix_rank(m, tol: nil, driver: 'svd') ⇒ Object

Compute matrix rank of array using SVD Rank is the number of singular values greater than tol.

Parameters:

  • m (Numo::NArray)

    matrix (>= 2-dimensional NArray)

  • tol (Float)

    threshold below which singular values are considered to be zero. If tol is nil, tol = sing_vals.max() * m.shape.max * EPSILON.

  • driver (String or Symbol)

    choose LAPACK solver from ‘svd’, ‘sdd’. (optional, default=’svd’)



765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
# File 'lib/numo/linalg/function.rb', line 765

def matrix_rank(m, tol:nil, driver:'svd')
  m = Numo::NArray.asarray(m)
  if m.ndim < 2
    m.ne(0).any? ? 1 : 0
  else
    case driver.to_s
    when /^(ge)?sdd$/, "turbo"
      s = Lapack.call(:gesdd, m, jobz:'N')[0]
    when /^(ge)?svd$/
      s = Lapack.call(:gesvd, m, jobu:'N', jobvt:'N')[0]
    else
      raise ArgumentError, "invalid driver: #{driver}"
    end
    tol ||= s.max(axis:-1, keepdims:true) *
      (m.shape[-2..-1].max * s.class::EPSILON)
    (s > tol).count(axis:-1)
  end
end

.norm(a, ord = nil, axis: nil, keepdims: false) ⇒ Numo::NArray

Compute matrix or vector norm.

|  ord  |  matrix norm           | vector norm                 |
| ----- | ---------------------- | --------------------------- |
|  nil  | Frobenius norm         | 2-norm                      |
| 'fro' | Frobenius 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) |

Parameters:

  • a (Numo::NArray)

    matrix or vector (>= 1-dimensinal NArray)

  • ord (String or Symbol) (defaults to: nil)

    Order of the norm .

  • axis (Integer or Array)

    Applied axes (optional).

  • keepdims (Bool)

    If true, the applied axes are left in result with size one (optional).

Returns:

  • (Numo::NArray)

    norm result



563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
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
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
# File 'lib/numo/linalg/function.rb', line 563

def norm(a, ord=nil, axis:nil, keepdims:false)
  a = Numo::NArray.asarray(a)

  # check axis
  if axis
    case axis
    when Integer
      axis = [axis]
    when Array
      if axis.size < 1 || axis.size > 2
        raise ArgmentError, "axis option should be 1- or 2-element array"
      end
    else
      raise ArgumentError, "invalid option for axis: #{axis}"
    end
    # swap axes
    if a.ndim > 1
      idx = (0...a.ndim).to_a
      tmp = []
      axis.each do |i|
        x = idx[i]
        if x.nil?
          raise ArgmentError, "axis contains same dimension"
        end
        tmp << x
        idx[i] = nil
      end
      idx.compact!
      idx.concat(tmp)
      a = a.transpose(*idx)
    end
  else
    case a.ndim
    when 0
      raise ArgumentError, "zero-dimensional array"
    when 1
      axis = [-1]
    else
      axis = [-2,-1]
    end
  end

  # calculate norm
  case axis.size

  when 1  # vector
    k = keepdims
    ord ||= 2  # default
    case ord.to_s
    when "0"
      r = a.class.cast(a.ne(0)).sum(axis:-1, keepdims:k)
    when "1"
      r = a.abs.sum(axis:-1, keepdims:k)
    when "2"
      r = Blas.call(:nrm2, a, keepdims:k)
    when /^-?\d+$/
      o = ord.to_i
      r = (a.abs**o).sum(axis:-1, keepdims:k)**(1.0/o)
    when /^inf(inity)?$/i
      r = a.abs.max(axis:-1, keepdims:k)
    when /^-inf(inity)?$/i
      r = a.abs.min(axis:-1, keepdims:k)
    else
      raise ArgumentError, "ord (#{ord}) is invalid for vector norm"
    end

  when 2  # matrix
    if keepdims
      fixdims = [true] * a.ndim
      axis.each do |i|
        if i < -a.ndim || i >= a.ndim
          raise ArgmentError, "axis (%d) is out of range", i
        end
        fixdims[i] = :new
      end
    end
    ord ||= "fro"  # default
    case ord.to_s
    when "1"
      r, = Lapack.call(:lange, a, '1')
    when "-1"
      r = a.abs.sum(axis:-2).min(axis:-1)
    when "2"
      svd, = Lapack.call(:gesvd, a, jobu:'N', jobvt:'N')
      r = svd.max(axis:-1)
    when "-2"
      svd, = Lapack.call(:gesvd, a, jobu:'N', jobvt:'N')
      r = svd.min(axis:-1)
    when /^f(ro)?$/i
      r, = Lapack.call(:lange, a, 'F')
    when /^inf(inity)?$/i
      r, = Lapack.call(:lange, a, 'I')
    when /^-inf(inity)?$/i
      r = a.abs.sum(axis:-1).min(axis:-1)
    else
      raise ArgumentError, "ord (#{ord}) is invalid for matrix norm"
    end
    if keepdims
      if NArray===r
        r = r[*fixdims]
      else
        r = a.class.new(1,1).store(r)
      end
    end
  end
  return r
end

.pinv(a, driver: "svd", rcond: nil) ⇒ Numo::NArray

Compute the (Moore-Penrose) pseudo-inverse of a matrix using svd or lstsq.

Examples:

a = Numo::DFloat.new(5,3).rand_norm
=> Numo::DFloat#shape=[5,3]
[[-0.581255, -0.168354, 0.586895],
 [-0.595142, -0.802802, -0.326106],
 [0.282922, 1.68427, 0.918499],
 [-0.0485384, -0.464453, -0.992194],
 [0.413794, -0.60717, -0.699695]]
b = Numo::Linalg.pinv(a,driver:"svd")
=> Numo::DFloat(view)#shape=[3,5]
[[-0.360863, -0.813125, -0.353367, -0.891963, 0.877253],
 [-0.227645, 0.162939, 0.696655, 0.787685, -0.469346],
 [0.408671, -0.308323, -0.337807, -1.13833, 0.228051]]
(a-a.dot(b.dot(a))).abs.max
=> 5.551115123125783e-16

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK driver from SVD (‘svd’, ‘sdd’) or Least square (‘lsd’,’lss’,’lsy’) (optional, default=’svd’)

  • rcond (Float)

    (optional, default=-1) RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

Returns:

  • (Numo::NArray)


963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
# File 'lib/numo/linalg/function.rb', line 963

def pinv(a, driver:"svd", rcond:nil)
  a = NArray.asarray(a)
  if a.ndim < 2
    raise NArray::ShapeError, "2-d array is required"
  end
  case driver
  when /^(ge)?s[dv]d$/
    s, u, vh = svd(a, driver:driver, job:'S')
    if rcond.nil? || rcond < 0
      rcond = ((SFloat===s) ? 1e3 : 1e6) * s.class::EPSILON
    elsif ! Numeric === rcond
      raise ArgumentError, "rcond must be Numeric"
    end
    cond = (s > rcond * s.max(axis:-1, keepdims:true))
    if cond.all?
      r = s.reciprocal
    else
      r = s.new_zeros
      r[cond] = s[cond].reciprocal
    end
    u *= r[false,:new,true]
    dot(u,vh).conj.swapaxes(-2,-1)
  when /^(ge)?ls[dsy]$/
    b = a.class.eye(a.shape[-2])
    x, = lstsq(a, b, driver:driver, rcond:rcond)
    x
  else
    raise ArgumentError, "#{driver.inspect} is not one of drivers: "+
      "svd, sdd, lsd, lss, lsy"
  end
end

.qr(a, mode: "reduce") ⇒ r, ...

Computes a QR factorization of a complex M-by-N matrix A: A = Q * R.

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • mode (String)
    • “reduce” – returns both Q and R,
    • “r” – returns only R,
    • “economy” – returns both Q and R but computed in economy-size,
    • “raw” – returns QR and TAU used in LAPACK.

Returns:

  • (r)

    if mode:”r”

  • ([q,r])

    if mode:”reduce” or “economic”

  • ([qr,tau])

    if mode:”raw” (LAPACK geqrf result)



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
# File 'lib/numo/linalg/function.rb', line 205

def qr(a, mode:"reduce")
  qr,tau, = Lapack.call(:geqrf, a)
  *shp,m,n = qr.shape
  r = (m >= n && %w[economic raw].include?(mode)) ?
    qr[false, 0...n, true].triu : qr.triu
  mode = mode.to_s.downcase
  case mode
  when "r"
    return r
  when "raw"
    return [qr,tau]
  when "reduce","economic"
    # skip
  else
    raise ArgumentError, "invalid mode:#{mode}"
  end
  if m < n
    q, = Lapack.call(:orgqr, qr[false, 0...m], tau)
  elsif mode == "economic"
    q, = Lapack.call(:orgqr, qr, tau)
  else
    qqr = qr.class.zeros(*(shp+[m,m]))
    qqr[false,0...n] = qr
    q, = Lapack.call(:orgqr, qqr, tau)
  end
  return [q,r]
end

.slogdet(a) ⇒ [sign,logdet]

Natural logarithm of the determinant of a matrix

Parameters:

  • a (Numo::NArray)

    matrix (>= 2-dimensional NArray)

Returns:

  • ([sign,logdet])
    • sign – A number representing the sign of the determinant.
    • logdet – The natural log of the absolute value of the determinant.


740
741
742
743
744
745
746
747
748
749
750
751
752
753
# File 'lib/numo/linalg/function.rb', line 740

def slogdet(a)
  lu, piv, = Lapack.call(:getrf, a)
  idx = piv.new_narray.store(piv.class.new(piv.shape[-1]).seq(1))
  m = piv.eq(idx).count_false(axis:-1) % 2
  sign = m * -2 + 1

  lud = lu.diagonal
  if (lud.eq 0).any?
    return 0, (-Float::INFINITY)
  end
  lud_abs = lud.abs
  sign *= (lud/lud_abs).prod
  [sign, NMath.log(lud_abs).sum(axis:-1)]
end

.solve(a, b, driver: "gen", uplo: 'U') ⇒ Numo::NArray

Solves linear equation a * x = b for x from square matrix a

Parameters:

  • a (Numo::NArray)

    n-by-n square matrix (>= 2-dimensinal NArray)

  • b (Numo::NArray)

    n-by-nrhs right-hand-side matrix (>= 1-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK diriver from ‘gen’,’sym’,’her’ or ‘pos’. (optional, default=’gen’)

  • uplo (String or Symbol)

    optional, default=’U’. Access upper or (‘U’) lower (‘L’) triangle. (omitted when driver:”gen”)

Returns:

  • (Numo::NArray)

    The solusion matrix/vector X.



798
799
800
801
802
803
804
805
806
807
808
809
# File 'lib/numo/linalg/function.rb', line 798

def solve(a, b, driver:"gen", uplo:'U')
  case driver.to_s
  when /^gen?(sv)?$/i
    # returns lu, x, ipiv, info
    Lapack.call(:gesv, a, b)[1]
  when /^(sym?|her?|pos?)(sv)?$/i
    func = driver[0..2].downcase+"sv"
    Lapack.call(func, a, b, uplo:uplo)[1]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.svd(a, driver: 'svd', job: 'A') ⇒ [sigma,u,vt]

Computes the Singular Value Decomposition (SVD) of a M-by-N matrix A, and the left and/or right singular vectors. The SVD is written

A = U * SIGMA * transpose(V)

where SIGMA is an M-by-N matrix which is zero except for its min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and V is an N-by-N orthogonal matrix. The diagonal elements of SIGMA are the singular values of A; they are real and non-negative, and are returned in descending order. The first min(m,n) columns of U and V are the left and right singular vectors of A. Note that the routine returns V**T, not V.

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK solver from ‘svd’, ‘sdd’. (optional, default=’svd’)

  • job (String or Symbol)
    • ‘A’: all M columns of U and all N rows of V**T are returned in the arrays U and VT.
    • ‘S’: the first min(M,N) columns of U and the first min(M,N) rows of V**T are returned in the arrays U and VT.
    • ‘N’: no columns of U or rows of V**T are computed.

Returns:

  • ([sigma,u,vt])

    SVD result. Array



258
259
260
261
262
263
264
265
266
267
268
269
270
# File 'lib/numo/linalg/function.rb', line 258

def svd(a, driver:'svd', job:'A')
  unless /^[ASN]/i =~ job
    raise ArgumentError, "invalid job: #{job.inspect}"
  end
  case driver.to_s
  when /^(ge)?sdd$/i, "turbo"
    Lapack.call(:gesdd, a, jobz:job)[0..2]
  when /^(ge)?svd$/i
    Lapack.call(:gesvd, a, jobu:job, jobvt:job)[0..2]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end

.svdvals(a, driver: 'svd') ⇒ Numo::NArray

Computes the Singular Values of a M-by-N matrix A. The SVD is written

A = U * SIGMA * transpose(V)

where SIGMA is an M-by-N matrix which is zero except for its min(m,n) diagonal elements. The diagonal elements of SIGMA are the singular values of A; they are real and non-negative, and are returned in descending order.

Parameters:

  • a (Numo::NArray)

    m-by-n matrix A (>= 2-dimensinal NArray)

  • driver (String or Symbol)

    choose LAPACK solver from ‘svd’, ‘sdd’. (optional, default=’svd’)

Returns:

  • (Numo::NArray)

    returns SIGMA (singular values).



287
288
289
290
291
292
293
294
295
296
# File 'lib/numo/linalg/function.rb', line 287

def svdvals(a, driver:'svd')
  case driver.to_s
  when /^(ge)?sdd$/i, "turbo"
    Lapack.call(:gesdd, a, jobz:'N')[0]
  when /^(ge)?svd$/i
    Lapack.call(:gesvd, a, jobu:'N', jobvt:'N')[0]
  else
    raise ArgumentError, "invalid driver: #{driver}"
  end
end