Fork me on GitHub

src/arraymancer/linear_algebra/helpers/decomposition_lapack

  Source Edit

Types

SupportedDecomposition = SomeFloat | Complex64 | Complex32
  Source Edit

Procs

proc geqrf[T: SupportedDecomposition](Q: var Tensor[T]; tau: var seq[T];
                                      scratchspace: var seq[T])

Wrapper for LAPACK geqrf routine (GEneral QR Factorization) Decomposition is done through Householder Reflection and without pivoting

In-place version, this will overwrite Q and tau

  Source Edit
proc gesdd[T: SupportedDecomposition; X: SupportedDecomposition](
    a: var Tensor[T]; U: var Tensor[T]; S: var Tensor[X]; Vh: var Tensor[T];
    scratchspace: var seq[T])

Wrapper for LAPACK gesdd routine (GEneral Singular value Decomposition by Divide & conquer)

Parameters:

  • a: Input - MxN matrix to factorize, in column major format
  • U: Output - Unitary matrix containing the left singular vectors as columns
  • S: Output - Singular values sorted in decreasing order
  • Vh: Output - Unitary matrix containing the right singular vectors as rows

SVD solves the equation: A = U S V.h

  • with S being a diagonal matrix of singular values
  • with V being the right singular vectors and V.h being the hermitian conjugate of V for real matrices, this is equivalent to V.t (transpose)

⚠️: Input must not contain NaN

Performance note:

  • Lapack, especially with the OpenBLAS backend is much more optimized for input M, N where M > N versus N < M (2x - 3x speed difference) Transpose accordingly. Matrices must be column major.
  Source Edit
proc getrf[T: SupportedDecomposition](lu: var Tensor[T];
                                      pivot_indices: var seq[int32])

Wrapper for LAPACK getrf routine (GEneral ??? Pivoted LU Factorization)

In-place version, this will overwrite LU and tau

  Source Edit
proc syevr[T: SupportedDecomposition](a: var Tensor[T]; uplo: static char;
                                      return_eigenvectors: static bool;
                                      low_idx: int; high_idx: int;
                                      eigenval, eigenvec: var Tensor[T];
                                      scratchspace: var seq[T])

Wrapper for LAPACK syevr routine (Symmetric Recursive Eigenvalue Decomposition)

eigenvalues are returned in ascending order (from lower to upper)

if uplo = 'L', the lower part of A is used it is destroyed on exit (and upper part is untouched) vice-versa if uplo = 'U' for the upper part of A

  Source Edit

Templates

template address(x: typed): untyped
  Source Edit
Arraymancer Technical reference Tutorial Spellbook (How-To's) Under the hood