Skip to content

blas #

The V Basic Linear Algebra System

This package implements Basic Linear Algebra System (BLAS) routines in V.

Backend Description Status Compilation Flags
BLAS Pure V implementation - High performance, zero dependencies Stable NONE
OpenBLAS OpenBLAS is an optimized BLAS library based on https://github.com/xianyi/OpenBLAS. Check the section OpenBLAS Backend for more information. Stable -d vsl_blas_cblas

🎉 Pure V Implementation

The pure V BLAS implementation is stable and production-ready, offering:

  • ✅ Zero Dependencies: No external C libraries required
  • ✅ High Performance: Competitive performance with C backends
  • ✅ Cross-Platform: Works on all platforms supported by V
  • ✅ Complete Coverage: All essential BLAS Level 1, 2, and 3 operations

Performance

The pure V implementation delivers excellent performance across a wide range of problem sizes. Benchmark results demonstrate competitive performance with optimized C backends while maintaining zero-dependency operation.

Run benchmarks to see performance characteristics:

v run benchmarks/blas_bench.v

Available Functions

Level 1 BLAS:

  • ddot - Dot product
  • dnrm2 - Euclidean norm
  • dasum - Sum of absolute values
  • daxpy - Vector addition (y = alpha*x + y)
  • dscal - Vector scaling (x = alpha*x)
  • dcopy - Vector copy
  • dswap - Vector swap
  • idamax - Index of maximum absolute value

Level 2 BLAS:

  • dgemv - General matrix-vector multiply
  • dger - Rank-1 update
  • dsymv - Symmetric matrix-vector multiply
  • dsyr - Symmetric rank-1 update
  • dtrmv - Triangular matrix-vector multiply
  • dtrsv - Triangular solve

Level 3 BLAS:

  • dgemm - General matrix-matrix multiply
  • dsymm - Symmetric matrix-matrix multiply
  • dsyrk - Symmetric rank-k update
  • dsyr2k - Symmetric rank-2k update
  • dtrmm - Triangular matrix-matrix multiply
  • dtrsm - Triangular solve with multiple RHS

Therefore, its routines are a little more lower level than the ones in the package vsl.la.

OpenBLAS Backend

We provide a backend for the OpenBLAS library. This backend is probably the fastest one for all platforms but it requires the installation of the OpenBLAS library.

Use the compilation flag -d vsl_blas_cblas to use the OpenBLAS backend instead of the pure V implementation and make sure that the OpenBLAS library is installed in your system.

Check the section below for more information about installing the OpenBLAS library.

Install dependencies

Homebrew (macOS)

brew install openblas

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 blas-openblas.

sudo pacman -S blas-openblas

fn cm_daxpy #

fn cm_daxpy(n int, alpha f64, x []f64, incx int, mut y []f64, incy int)

cm_daxpy exposes this operation as part of the public API.

fn cm_ddot #

fn cm_ddot(n int, x []f64, incx int, y []f64, incy int) f64

cm_ddot exposes this operation as part of the public API.

fn cm_dgemm #

fn cm_dgemm(trans_a Transpose, trans_b Transpose, m int, n int, k int, alpha f64, a []f64, lda int, b []f64, ldb int, beta f64, mut c []f64, ldc int)

cm_dgemm exposes this operation as part of the public API.

fn cm_dgemv #

fn cm_dgemv(trans Transpose, m int, n int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)

cm_dgemv exposes this operation as part of the public API.

fn cm_dscal #

fn cm_dscal(n int, alpha f64, mut x []f64, incx int)

cm_dscal exposes this operation as part of the public API.

fn cm_dsyrk #

fn cm_dsyrk(uplo Uplo, trans Transpose, n int, k int, alpha f64, a []f64, lda int, beta f64, mut c []f64, ldc int)

cm_dsyrk exposes this operation as part of the public API.

fn cm_dtrsm #

fn cm_dtrsm(side Side, uplo Uplo, trans_a Transpose, diag Diagonal, m int, n int, alpha f64, a []f64, lda int, mut b []f64, ldb int)

cm_dtrsm exposes this operation as part of the public API.

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

dasum exposes this operation as part of the public API.

fn daxpy #

fn daxpy(n int, alpha f64, x []f64, incx int, mut y []f64, incy int)

daxpy exposes this operation as part of the public API.

fn dcopy #

fn dcopy(n int, x []f64, incx int, mut y []f64, incy int)

dcopy exposes this operation as part of the public API.

fn ddot #

fn ddot(n int, x []f64, incx int, y []f64, incy int) f64

ddot exposes this operation as part of the public API.

fn dgbmv #

fn dgbmv(trans_a Transpose, m int, n int, kl int, ku int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)

dgbmv exposes this operation as part of the public API.

fn dgemm #

fn dgemm(trans_a Transpose, trans_b Transpose, m int, n int, k int, alpha f64, a []f64, lda int, b []f64, ldb int, beta f64, mut cc []f64, ldc int)

dgemm exposes this operation as part of the public API.

fn dgemv #

fn dgemv(trans Transpose, m int, n int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)

dgemv exposes this operation as part of the public API.

fn dger #

fn dger(m int, n int, alpha f64, x []f64, incx int, y []f64, incy int, mut a []f64, lda int)

dger exposes this operation as part of the public API.

fn dnrm2 #

fn dnrm2(n int, x []f64, incx int) f64

dnrm2 exposes this operation as part of the public API.

fn drot #

fn drot(n int, mut x []f64, incx int, mut y []f64, incy int, c f64, s f64)

drot exposes this operation as part of the public API.

fn dsbmv #

fn dsbmv(uplo Uplo, n int, k int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)

dsbmv exposes this operation as part of the public API.

fn dscal #

fn dscal(n int, alpha f64, mut x []f64, incx int)

dscal exposes this operation as part of the public API.

fn dspmv #

fn dspmv(uplo Uplo, n int, alpha f64, ap []f64, x []f64, incx int, beta f64, mut y []f64, incy int)

dspmv exposes this operation as part of the public API.

fn dspr #

fn dspr(uplo Uplo, n int, alpha f64, x []f64, incx int, mut ap []f64)

dspr exposes this operation as part of the public API.

fn dspr2 #

fn dspr2(uplo Uplo, n int, alpha f64, x []f64, incx int, y []f64, incy int, mut ap []f64)

dspr2 exposes this operation as part of the public API.

fn dswap #

fn dswap(n int, mut x []f64, incx int, mut y []f64, incy int)

dswap exposes this operation as part of the public API.

fn dsymv #

fn dsymv(uplo Uplo, n int, alpha f64, a []f64, lda int, x []f64, incx int, beta f64, mut y []f64, incy int)

dsymv exposes this operation as part of the public API.

fn dsyr #

fn dsyr(uplo Uplo, n int, alpha f64, x []f64, incx int, mut a []f64, lda int)

dsyr exposes this operation as part of the public API.

fn dsyr2 #

fn dsyr2(uplo Uplo, n int, alpha f64, x []f64, incx int, y []f64, incy int, mut a []f64, lda int)

dsyr2 exposes this operation as part of the public API.

fn dsyr2k #

fn dsyr2k(uplo Uplo, trans_a Transpose, n int, k int, alpha f64, a []f64, lda int, b []f64, ldb int, beta f64, mut c []f64, ldc int)

dsyr2k exposes this operation as part of the public API.

fn dsyrk #

fn dsyrk(uplo Uplo, trans_a Transpose, n int, k int, alpha f64, a []f64, lda int, beta f64, mut c []f64, ldc int)

dsyrk exposes this operation as part of the public API.

fn dtbmv #

fn dtbmv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, k int, a []f64, lda int, mut x []f64, incx int)

dtbmv exposes this operation as part of the public API.

fn dtbsv #

fn dtbsv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, k int, a []f64, lda int, mut x []f64, incx int)

dtbsv exposes this operation as part of the public API.

fn dtpmv #

fn dtpmv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, ap []f64, mut x []f64, incx int)

dtpmv exposes this operation as part of the public API.

fn dtpsv #

fn dtpsv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, ap []f64, mut x []f64, incx int)

dtpsv exposes this operation as part of the public API.

fn dtrmm #

fn dtrmm(side Side, uplo Uplo, trans Transpose, diag Diagonal, m int, n int, alpha f64, a []f64, lda int, mut b []f64, ldb int)

dtrmm performs triangular matrix multiplication. Input matrices are expected in row-major format (as used by la/ module and tests). blas64.dtrmm uses row-major access pattern but validates ldb >= m (inconsistent). We pass matrices directly but ensure ldb >= m for validation.

fn dtrmv #

fn dtrmv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, a []f64, lda int, mut x []f64, incx int)

dtrmv exposes this operation as part of the public API.

fn dtrsm #

fn dtrsm(side Side, uplo Uplo, trans Transpose, diag Diagonal, m int, n int, alpha f64, a []f64, lda int, mut b []f64, ldb int)

dtrsm exposes this operation as part of the public API.

fn dtrsv #

fn dtrsv(uplo Uplo, trans_a Transpose, diag Diagonal, n int, a []f64, lda int, mut x []f64, incx int)

dtrsv exposes this operation as part of the public API.

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 from_blas64_diagonal #

fn from_blas64_diagonal(diag blas64.Diagonal) Diagonal

Convert blas64.Diagonal to VSL Diagonal

fn from_blas64_layout #

fn from_blas64_layout(layout blas64.MemoryLayout) MemoryLayout

Convert blas64.MemoryLayout to VSL MemoryLayout

fn from_blas64_side #

fn from_blas64_side(side blas64.Side) Side

Convert blas64.Side to VSL Side

fn from_blas64_transpose #

fn from_blas64_transpose(trans blas64.Transpose) Transpose

Convert blas64.Transpose to VSL Transpose

fn from_blas64_uplo #

fn from_blas64_uplo(uplo blas64.Uplo) Uplo

Convert blas64.Uplo to VSL Uplo

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 idamax #

fn idamax(n int, x []f64, incx int) int

idamax exposes this operation as part of the public API.

fn join_complex #

fn join_complex(v_real []f64, v_imag []f64) []complex.Complex

join_complex joins real and imag parts of array

fn set_num_threads #

fn set_num_threads(n int)

set_num_threads exposes this operation as part of the public API.

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

fn to_blas64_diagonal #

fn to_blas64_diagonal(diag Diagonal) blas64.Diagonal

Convert VSL Diagonal to blas64.Diagonal

fn to_blas64_layout #

fn to_blas64_layout(layout MemoryLayout) blas64.MemoryLayout

Convert VSL MemoryLayout to blas64.MemoryLayout

fn to_blas64_side #

fn to_blas64_side(side Side) blas64.Side

Convert VSL Side to blas64.Side

fn to_blas64_transpose #

fn to_blas64_transpose(trans Transpose) blas64.Transpose

Convert VSL Transpose to blas64.Transpose

fn to_blas64_uplo #

fn to_blas64_uplo(uplo Uplo) blas64.Uplo

Convert VSL Uplo to blas64.Uplo

fn Diagonal.from #

fn Diagonal.from[W](input W) !Diagonal

fn MemoryLayout.from #

fn MemoryLayout.from[W](input W) !MemoryLayout

fn Side.from #

fn Side.from[W](input W) !Side

fn Transpose.from #

fn Transpose.from[W](input W) !Transpose

fn Uplo.from #

fn Uplo.from[W](input W) !Uplo

enum Diagonal #

enum Diagonal {
	non_unit = 131
	unit     = 132
}

Diagonal is used to specify whether the diagonal of a matrix is unit or non-unit.

enum MemoryLayout #

enum MemoryLayout {
	row_major = 101
	col_major = 102
}

MemoryLayout is used to specify the memory layout of a matrix.

enum Side #

enum Side {
	left  = 141
	right = 142
}

Side is used to specify whether a matrix is on the left or right side in a matrix-matrix multiplication.

enum Transpose #

enum Transpose {
	no_trans      = 111
	trans         = 112
	conj_trans    = 113
	conj_no_trans = 114
}

Transpose is used to specify the transposition of a matrix.

enum Uplo #

enum Uplo {
	upper = 121
	lower = 122
	all   = 99
}

Uplo is used to specify whether the upper or lower triangle of a matrix is