vlas #
V Linear Algebra System
This package implements BLAS and LAPACKE functions. It provides different backends:
Backend | Description | Status | Compilation Flags |
---|---|---|---|
VLAS | Pure V implementation | WIP | NONE |
OpenBLAS | OpenBLAS is an optimized BLAS library based on https://github.com/xianyi/OpenBLAS. Check the section OpenBLAS Backend for more information. | Working | -d vsl_vlas_cblas |
LAPACKE | LAPACKE is a C interface to the LAPACK linear algebra routines | Working | -d vsl_vlas_lapacke |
Therefore, its routines are a little more lower level than the ones in the package vsl.la
.
OpenBLAS Backend
Use the flag -d vsl_vlas_cblas
to use the OpenBLAS backend.
Install dependencies
Debian/Ubuntu GNU Linux
libopenblas-dev
is not needed when using the pure V backend.
sudo apt-get install -y --no-install-recommends
Arch Linux/Manjaro GNU Linux
The best way of installing OpenBlas is using lapack-openblas.
yay -S lapack-openblas
or
git clone https://aur.archlinux.org/lapack-openblas.git /tmp/lapack-openblas
cd /tmp/lapack-openblas
makepkg -si
macOS
brew install openblas
LAPACKE Backend
Use the flag -d vsl_vlas_lapacke
to use the LAPACKE backend (enabled by default for now).
Install dependencies
Debian/Ubuntu GNU Linux
sudo apt-get install -y --no-install-recommends
Arch Linux/Manjaro GNU Linux
The best way of installing LAPACKE is using lapack-openblas.
yay -S lapack-openblas
or
git clone https://aur.archlinux.org/lapack-openblas.git /tmp/lapack-openblas
cd /tmp/lapack-openblas
makepkg -si
fn c_trans #
fn c_trans(trans bool) blas.Transpose
fn c_uplo #
fn c_uplo(up bool) blas.Uplo
fn col_major_complex_to_slice #
fn col_major_complex_to_slice(m int, n int, data []complex.Complex) [][]complex.Complex
col_major_complex_to_slice converts col-major matrix to nested slice
fn col_major_to_slice #
fn col_major_to_slice(m int, n int, data []f64) [][]f64
col_major_to_slice converts col-major matrix to nested slice
fn dasum #
fn dasum(n int, x []f64, incx int) f64
fn daxpy #
fn daxpy(n int, alpha f64, x []f64, incx int, mut y []f64, incy int)
fn dcopy #
fn dcopy(n int, x []f64, incx int, mut y []f64, incy int)
fn ddot #
fn ddot(n int, x []f64, incx int, y []f64, incy int) f64
fn dgeev #
fn dgeev(calc_vl bool, calc_vr bool, n int, mut a []f64, lda int, wr []f64, wi []f64, vl []f64, ldvl_ int, vr []f64, ldvr_ int)
dgeev computes for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
See: http://www.netlib.org/lapack/explore-html/d9/d28/dgeev_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-geev
See: https://www.nag.co.uk/numeric/fl/nagdoc_fl26/html/f08/f08naf.html
The right eigenvector v(j) of A satisfies
A * v(j) = lambda(j) * v(j)
where lambda(j) is its eigenvalue.
The left eigenvector u(j) of A satisfies
u(j)**H * A = lambda(j) * u(j)**H
where u(j)**H denotes the conjugate-transpose of u(j).
The computed eigenvectors are normalized to have Euclidean norm equal to 1 and largest component real.
fn dgemm #
fn dgemm(trans_a bool, trans_b bool, m int, n int, k int, alpha f64, a []f64, lda int, b []f64, ldb int, beta f64, mut cc []f64, ldc int)
fn dgemv #
fn dgemv(trans bool, m int, n int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)
fn dger #
fn dger(m int, n int, alpha f64, x []f64, incx int, y []f64, incy int, mut a []f64, lda int)
fn dgesv #
fn dgesv(n int, nrhs int, mut a []f64, lda int, ipiv []int, mut b []f64, ldb int)
dgesv computes the solution to a real system of linear equations.
See: http://www.netlib.org/lapack/explore-html/d8/d72/dgesv_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-gesv
The system is:
A * X = B,
where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
A = P * L * U,
where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.
Note: matrix 'a' will be modified
fn dgesvd #
fn dgesvd(jobu rune, jobvt rune, m int, n int, a []f64, lda int, s []f64, u []f64, ldu int, vt []f64, ldvt int, superb []f64)
dgesvd computes the singular value decomposition (SVD) of a real M-by-N matrix A, optionally computing the left and/or right singular vectors.
See: http://www.netlib.org/lapack/explore-html/d8/d2d/dgesvd_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-gesvd
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.
Note: matrix 'a' will be modified
fn dgetrf #
fn dgetrf(m int, n int, mut a []f64, lda int, ipiv []int)
dgetrf computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.
See: http://www.netlib.org/lapack/explore-html/d3/d6a/dgetrf_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-getrf
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).
This is the right-looking Level 3 BLAS version of the algorithm.
Note: (1) matrix 'a' will be modified(2) ipiv indices are 1-based (i.e. Fortran)
fn dgetri #
fn dgetri(n int, mut a []f64, lda int, ipiv []int)
dgetri computes the inverse of a matrix using the LU factorization computed by DGETRF.
See: http://www.netlib.org/lapack/explore-html/df/da4/dgetri_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-getri
This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).
fn dlange #
fn dlange(norm rune, m int, n int, a []f64, lda int, work []f64) f64
fn dnrm2 #
fn dnrm2(n int, x []f64, incx int) f64
fn dpotrf #
fn dpotrf(up bool, n int, mut a []f64, lda int)
dpotrf computes the Cholesky factorization of a real symmetric positive definite matrix A.
See: http://www.netlib.org/lapack/explore-html/d0/d8a/dpotrf_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-potrf
The factorization has the form
A = U**T * U, if UPLO = 'U'
or
A = L * L**T, if UPLO = 'L'
where U is an upper triangular matrix and L is lower triangular.
This is the block version of the algorithm, calling Level 3 BLAS.
fn drot #
fn drot(n int, mut x []f64, incx int, mut y []f64, incy int, c f64, s f64)
fn dscal #
fn dscal(n int, alpha f64, mut x []f64, incx int)
fn dswap #
fn dswap(n int, mut x []f64, incx int, mut y []f64, incy int)
fn dsyr #
fn dsyr(uplo bool, n int, alpha f64, x []f64, incx int, mut a []f64, lda int)
fn dsyr2 #
fn dsyr2(uplo bool, n int, alpha f64, x []f64, incx int, y []f64, incy int, mut a []f64, lda int)
fn dtrmv #
fn dtrmv(uplo bool, trans_a bool, diag blas.Diagonal, n int, a []f64, lda int, mut x []f64, incx int)
fn dtrsv #
fn dtrsv(uplo bool, trans_a bool, diag blas.Diagonal, n int, a []f64, lda int, mut x []f64, incx int)
fn eigenvecs_build #
fn eigenvecs_build(mut vv []complex.Complex, wr []f64, wi []f64, v []f64)
eigenvecs_build builds complex eigenvectros created by Dgeev function
input: wr
, wi
: real and imag parts of eigenvalues. v
: left or right eigenvectors from Dgeev.
output: vv
: complex version of left or right eigenvector [pre-allocated].
NOTE: (no checks made).
n = wr.len = wi.len = v.len
2 * n = vv.len
fn eigenvecs_build_both #
fn eigenvecs_build_both(mut vvl []complex.Complex, mut vvr []complex.Complex, wr []f64, wi []f64, vl []f64, vr []f64)
eigenvecs_build_both builds complex left and right eigenvectros created by Dgeev function
input: wr
, wi
:real and imag parts of eigenvalues. vl
, vr
:left and right eigenvectors from Dgeev.
output: vvl
, vvr
:complex version of left and right eigenvectors [pre-allocated].
NOTE: (no checks made).
n = wr.len = wi.len = vl.len = vr.len
2 * n = vvl.len = vvr.len
fn extract_col #
fn extract_col(j int, m int, n int, a []f64) []f64
extract_col extracts j column from (m,n) col-major matrix
fn extract_col_complex #
fn extract_col_complex(j int, m int, n int, a []complex.Complex) []complex.Complex
extract_col_complex extracts j column from (m,n) col-major matrix (complex version)
fn extract_row #
fn extract_row(i int, m int, n int, a []f64) []f64
extract_row extracts i row from (m,n) col-major matrix
fn extract_row_complex #
fn extract_row_complex(i int, m int, n int, a []complex.Complex) []complex.Complex
extract_row_complex extracts i row from (m,n) col-major matrix (complex version)
fn get_join_complex #
fn get_join_complex(v_real []f64, v_imag []f64) []complex.Complex
get_join_complex joins real and imag parts of array
fn get_split_complex #
fn get_split_complex(v []complex.Complex) ([]f64, []f64)
get_split_complex splits real and imag parts of array
fn join_complex #
fn join_complex(v_real []f64, v_imag []f64) []complex.Complex
join_complex joins real and imag parts of array
fn print_col_major #
fn print_col_major(m int, n int, data []f64, nfmt_ string) string
print_col_major prints matrix (without commas or brackets)
fn print_col_major_complex #
fn print_col_major_complex(m int, n int, data []complex.Complex, nfmt_r_ string, nfmt_i_ string) string
print_col_major_complex prints matrix (without commas or brackets). NOTE: if non-empty, nfmt_i must have '+' e.g. %+g
fn print_col_major_complex_v #
fn print_col_major_complex_v(m int, n int, data []complex.Complex, nfmt_r_ string, nfmt_i_ string) string
print_col_major_complex_v prints matrix in v format NOTE: if non-empty, nfmt_i must have '+' e.g. %+g
fn print_col_major_omplex_py #
fn print_col_major_omplex_py(m int, n int, data []complex.Complex, nfmt_r_ string, nfmt_i_ string) string
print_col_major_omplex_py prints matrix in Python format NOTE: if non-empty, nfmt_i must have '+' e.g. %+g
fn print_col_major_py #
fn print_col_major_py(m int, n int, data []f64, nfmt_ string) string
print_col_major_py prints matrix in Python format
fn print_col_major_v #
fn print_col_major_v(m int, n int, data []f64, nfmt_ string) string
print_col_major_v prints matrix in v format
fn set_num_threads #
fn set_num_threads(n int)
set_num_threads sets the number of threads in VLAS
fn slice_to_col_major #
fn slice_to_col_major(a [][]f64) []f64
slice_to_col_major converts nested slice into an array representing a col-major matrix
NOTE: make sure to have at least 1x1 item
fn slice_to_col_major_complex #
fn slice_to_col_major_complex(a [][]complex.Complex) []complex.Complex
slice_to_col_major_complex converts nested slice into an array representing a col-major matrix of complex numbers.
data[i+j*m] = a[i][j]
NOTE: make sure to have at least 1x1 item
fn split_complex #
fn split_complex(v []complex.Complex) ([]f64, []f64)
split_complex splits real and imag parts of array
- README
- fn c_trans
- fn c_uplo
- fn col_major_complex_to_slice
- fn col_major_to_slice
- fn dasum
- fn daxpy
- fn dcopy
- fn ddot
- fn dgeev
- fn dgemm
- fn dgemv
- fn dger
- fn dgesv
- fn dgesvd
- fn dgetrf
- fn dgetri
- fn dlange
- fn dnrm2
- fn dpotrf
- fn drot
- fn dscal
- fn dswap
- fn dsyr
- fn dsyr2
- fn dtrmv
- fn dtrsv
- fn eigenvecs_build
- fn eigenvecs_build_both
- fn extract_col
- fn extract_col_complex
- fn extract_row
- fn extract_row_complex
- fn get_join_complex
- fn get_split_complex
- fn join_complex
- fn print_col_major
- fn print_col_major_complex
- fn print_col_major_complex_v
- fn print_col_major_omplex_py
- fn print_col_major_py
- fn print_col_major_v
- fn set_num_threads
- fn slice_to_col_major
- fn slice_to_col_major_complex
- fn split_complex