The stdlib_sparse_spmv submodule provides matrix-vector product kernels.
#:include "common.fypp" #:set RANKS = range(1, 2+1) #:set R_KINDS_TYPES = list(zip(REAL_KINDS, REAL_TYPES, REAL_SUFFIX)) #:set C_KINDS_TYPES = list(zip(CMPLX_KINDS, CMPLX_TYPES, CMPLX_SUFFIX)) #:set KINDS_TYPES = R_KINDS_TYPES+C_KINDS_TYPES #! define ranks without parentheses #:def rksfx2(rank) #{if rank > 0}#${":," + ":," * (rank - 1)}$#{endif}# #:enddef !! The `stdlib_sparse_spmv` submodule provides matrix-vector product kernels. !! ! This code was modified from https://github.com/jalvesz/FSPARSE by its author: Alves Jose module stdlib_sparse_spmv use stdlib_sparse_constants use stdlib_sparse_kinds implicit none private !! Version experimental !! !! Apply the sparse matrix-vector product $$y = \alpha * op(M) * x + \beta * y $$ !! [Specifications](../page/specs/stdlib_sparse.html#spmv) interface spmv #:for k1, t1, s1 in (KINDS_TYPES) #:for rank in RANKS module subroutine spmv_coo_${rank}$d_${s1}$(matrix,vec_x,vec_y,alpha,beta,op) type(COO_${s1}$_type), intent(in) :: matrix ${t1}$, intent(in) :: vec_x${ranksuffix(rank)}$ ${t1}$, intent(inout) :: vec_y${ranksuffix(rank)}$ ${t1}$, intent(in), optional :: alpha ${t1}$, intent(in), optional :: beta character(1), intent(in), optional :: op end subroutine module subroutine spmv_csr_${rank}$d_${s1}$(matrix,vec_x,vec_y,alpha,beta,op) type(CSR_${s1}$_type), intent(in) :: matrix ${t1}$, intent(in) :: vec_x${ranksuffix(rank)}$ ${t1}$, intent(inout) :: vec_y${ranksuffix(rank)}$ ${t1}$, intent(in), optional :: alpha ${t1}$, intent(in), optional :: beta character(1), intent(in), optional :: op end subroutine module subroutine spmv_csc_${rank}$d_${s1}$(matrix,vec_x,vec_y,alpha,beta,op) type(CSC_${s1}$_type), intent(in) :: matrix ${t1}$, intent(in) :: vec_x${ranksuffix(rank)}$ ${t1}$, intent(inout) :: vec_y${ranksuffix(rank)}$ ${t1}$, intent(in), optional :: alpha ${t1}$, intent(in), optional :: beta character(1), intent(in), optional :: op end subroutine module subroutine spmv_ell_${rank}$d_${s1}$(matrix,vec_x,vec_y,alpha,beta,op) type(ELL_${s1}$_type), intent(in) :: matrix ${t1}$, intent(in) :: vec_x${ranksuffix(rank)}$ ${t1}$, intent(inout) :: vec_y${ranksuffix(rank)}$ ${t1}$, intent(in), optional :: alpha ${t1}$, intent(in), optional :: beta character(1), intent(in), optional :: op end subroutine #:endfor module subroutine spmv_sellc_${s1}$(matrix,vec_x,vec_y,alpha,beta,op) type(SELLC_${s1}$_type), intent(in) :: matrix ${t1}$, intent(in) :: vec_x(:) ${t1}$, intent(inout) :: vec_y(:) ${t1}$, intent(in), optional :: alpha ${t1}$, intent(in), optional :: beta character(1), intent(in), optional :: op end subroutine #:endfor end interface public :: spmv end module