created fortran interface for C spmm code

sp3mm-interface
wlthr 2 years ago
parent 0c88352530
commit 250163f1bc

@ -1,44 +1,54 @@
module sp3mm_mod
use iso_c_binding
use psb_const_mod
use psb_error_mod
! use psb_const_mod
! use psb_error_mod
interface spmm_row_by_row
subroutine dspmm_row_by_row_ub(a,b,c,info)
use psb_d_mat_mod, only : psb_dspmat_type
interface psb_dspmm
subroutine dspmm(a,b,c,info)
use psb_mat_mod
import :: psb_ipk_
implicit none
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out) :: c
integer(psb_ipk_), intent(out) :: info
end subroutine dspmm_row_by_row_ub
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out):: c
integer(psb_ipk_), intent(out) :: info
end subroutine dspmm
end interface psb_dspmm
subroutine dspmm_row_by_row_symb_num(a,b,c,info)
use psb_d_mat_mod, only : psb_dspmat_type
import :: psb_ipk_
implicit none
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out) :: c
integer(psb_ipk_), intent(out) :: info
end subroutine dspmm_row_by_row_symb_num
! interface spmm_row_by_row
! subroutine dspmm_row_by_row_ub(a,b,c,info)
! use psb_d_mat_mod, only : psb_dspmat_type
! import :: psb_ipk_
! implicit none
! type(psb_dspmat_type), intent(in) :: a,b
! type(psb_dspmat_type), intent(out) :: c
! integer(psb_ipk_), intent(out) :: info
! end subroutine dspmm_row_by_row_ub
subroutine dspmm_row_by_row_1d_blocks_symb_num(a,b,c,info)
use psb_d_mat_mod, only : psb_dspmat_type
import :: psb_ipk_
implicit none
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out) :: c
integer(psb_ipk_), intent(out) :: info
end subroutine dspmm_row_by_row_1d_blocks_symb_num
! subroutine dspmm_row_by_row_symb_num(a,b,c,info)
! use psb_d_mat_mod, only : psb_dspmat_type
! import :: psb_ipk_
! implicit none
! type(psb_dspmat_type), intent(in) :: a,b
! type(psb_dspmat_type), intent(out) :: c
! integer(psb_ipk_), intent(out) :: info
! end subroutine dspmm_row_by_row_symb_num
subroutine dspmm_row_by_row_2d_blocks_symb_num(a,b,c,info)
use psb_d_mat_mod, only : psb_dspmat_type
import :: psb_ipk_
implicit none
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out) :: c
integer(psb_ipk_), intent(out) :: info
end subroutine dspmm_row_by_row_2d_blocks_symb_num
end interface spmm_row_by_row
! subroutine dspmm_row_by_row_1d_blocks_symb_num(a,b,c,info)
! use psb_d_mat_mod, only : psb_dspmat_type
! import :: psb_ipk_
! implicit none
! type(psb_dspmat_type), intent(in) :: a,b
! type(psb_dspmat_type), intent(out) :: c
! integer(psb_ipk_), intent(out) :: info
! end subroutine dspmm_row_by_row_1d_blocks_symb_num
! subroutine dspmm_row_by_row_2d_blocks_symb_num(a,b,c,info)
! use psb_d_mat_mod, only : psb_dspmat_type
! import :: psb_ipk_
! implicit none
! type(psb_dspmat_type), intent(in) :: a,b
! type(psb_dspmat_type), intent(out) :: c
! integer(psb_ipk_), intent(out) :: info
! end subroutine dspmm_row_by_row_2d_blocks_symb_num
! end interface spmm_row_by_row
end module sp3mm_mod

@ -1,6 +1,8 @@
#include "../include/Sp3MM_CSR_OMP_Multi.h"
#include "../include/utils.h"
enum impl_types {ROW_BY_ROW_UB};
/**
* @brief performs multiplication of two sparse matrices
* A and B stored in the CRS format. The resulting matrix is C
@ -25,40 +27,43 @@
* each row in b_as and b_ja
* @param[in] b_rl array with the lengths of each rows of B
* @param[in] b_max_row_nz maximum number of coefficients in a row in B
* @param[out] c_m number of colums of C
* @param[out] c_n number of rows of C
* @param[out] c_nz number of non zero elements in C
* @param[out] c_as array with the non zero coefficients of C
* @param[out] c_ja array with the column indices of the non
* zero coefficients of C
* @param[out] c_irp array with the indices of the beginning of
* each row in c_as and c_ja
* @param[out] c_rl array with the lengths of each rows of C
* @param[out] c_max_row_nz maximum number of coefficients in a row in C
* @param[out] info return value to check if the operation was successful
* @param[in] impl_choice implementation choice
* @param[out] c_m number of colums of C
* @param[out] c_n number of rows of C
* @param[out] c_nz number of non zero elements in C
* @param[out] c_as array with the non zero coefficients of C
* @param[out] c_ja array with the column indices of the non
* zero coefficients of C
* @param[out] c_irp array with the indices of the beginning of
* each row in c_as and c_ja
* @param[out] c_rl array with the lengths of each rows of C
* @param[out] c_max_row_nz maximum number of coefficients in a row in C
* @param[out] info return value to check if the operation was successful
*/
#ifdef ROWLENS
void psb_f_spmm_row_by_row_ub_0(idx_t a_m, idx_t a_n, idx_t a_nz,
void psb_f_spmm(idx_t a_m, idx_t a_n, idx_t a_nz,
double *a_as, idx_t *a_ja,
idx_t *a_irp, idx_t *a_rl, idx_t a_max_row_nz,
idx_t b_m, idx_t b_n, idx_t b_nz,
double *b_as, idx_t *b_ja,
idx_t *b_irp, idx_t *b_rl, idx_t b_max_row_nz,
enum impl_types impl_choice,
idx_t *c_m, idx_t *c_n, idx_t *c_nz,
double **c_as, idx_t **c_ja,
idx_t **c_irp, idx_t **c_rl, idx_t *c_max_row_nz,
int info)
#else
void psb_f_spmm_row_by_row_ub_0(idx_t a_m, idx_t a_n, idx_t a_nz,
void psb_f_spmm(idx_t a_m, idx_t a_n, idx_t a_nz,
double *a_as, idx_t *a_ja,
idx_t *a_irp, idx_t a_max_row_nz,
idx_t b_m, idx_t b_n, idx_t b_nz,
double *b_as, idx_t *b_ja,
idx_t *b_irp, idx_t b_max_row_nz,
enum impl_types impl_choice,
idx_t *c_m, idx_t *c_n, idx_t *c_nz,
double **c_as, idx_t **c_ja,
idx_t **c_irp, idx_t *c_max_row_nz,
int info)
int *info)
#endif
{
int rc;
@ -92,7 +97,14 @@ void psb_f_spmm_row_by_row_ub_0(idx_t a_m, idx_t a_n, idx_t a_nz,
b->MAX_ROW_NZ = b_max_row_nz;
// performing spmm
c = spmmRowByRow_0(a, b, cfg);
switch (impl_choice)
{
case ROW_BY_ROW_UB:
c = spmmRowByRow_0(a, b, cfg);
break;
default:
break;
}
// output result
*(c_m) = c->M;

@ -1,12 +1,95 @@
subroutine dspmm_row_by_row_ub(a,b,c,info)
use psb_error_mod
use psb_base_mat_mod
use psb_d_mat_mod, only : psb_dspmat_type
use psb_objhandle_mod, only: spmat_t, config_t
implicit none
type(psb_dspmat_type), intent(in) :: a,b
type(psb_dspmat_type), intent(out) :: c
integer(psb_ipk_), intent(out) :: info
! TODO : implement the C interface
end subroutine dspmm_row_by_row_ub
subroutine dspmm(a,b,c,info, impl_choice)
use psb_d_mat_mod
implicit none
type(psb_d_csr_sparse_mat), intent(in) :: a,b
type(psb_d_csr_sparse_mat), intent(out):: c
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), intent(in), optional :: impl_choice
! Internal variables
integer(c_size_t), value :: a_m,a_n,a_nz
real(c_double) :: a_as
integer(c_size_t) :: a_ja,a_irp
integer(c_size_t), value :: a_max_row_nz
integer(c_size_t), value :: b_m,b_n,b_nz
real(c_double) :: b_as
integer(c_size_t) :: b_ja,b_irp
integer(c_size_t), value :: b_max_row_nz
integer(c_int), value :: impl_choice_
integer(c_size_t), value :: c_m,c_n,c_nz
real(c_double) :: c_as
integer(c_size_t) :: c_ja,c_irp
integer(c_size_t), value :: c_max_row_nz
! Initializing internal variables
a_m = a%get_nrows()
a_n = a%get_ncols()
a_nz = a%get_nzeros()
a_as = a%val
a_ja = a%ja
a_irp = a%irp
! a_max_row_nz
b_m = b%get_nrows()
b_n = b%get_ncols()
b_nz = b%get_nzeros()
b_as = b%val
b_ja = b%ja
b_irp = b%irp
if (present(impl_choice)) then
impl_choice_ = impl_choice
else
impl_choice_ = 0
end if
! Calling psb_f_spmm
psb_f_spmm(a_m,a_n,a_nz,&
& a_as,a_ja,a_irp,&
& a_max_row_nz,&
& b_m,b_n,b_nz,&
& b_as,b_ja,b_irp,&
& impl_choice_,&
& c_m,c_n,c_nz,&
& c_as,c_ja,c_irp,&
& c_max_row_nz,&
& info)
! Putting results in c
c%m = c_m
c%n = c_n
c%val = c_as
c%ja = c_ja
c%irp = c_irp
interface
subroutine psb_f_spmm(c_a_m,c_a_n,c_a_nz,&
& c_a_as,c_a_ja,c_a_irp,&
& c_a_max_row_nz,
& c_b_m,c_b_n,c_b_nz,&
& c_b_as,c_b_ja,c_b_irp,&
& c_b_max_row_nz,
& c_impl_choice,&
& c_c_m,c_c_n,c_c_nz,&
& c_c_as,c_c_ja,c_c_irp,
& c_c_max_row_nz,&
& c_info) bind(c)
use iso_c_binding
use psb_base_mod
integer(c_size_t), intent(in), value :: c_a_m,c_a_n,c_a_nz
real(c_double), intent(in) :: c_a_as
integer(c_size_t), intent(in) :: c_a_ja,c_a_irp
integer(c_size_t), intent(in), value :: c_a_max_row_nz
integer(c_size_t), intent(in), value :: c_b_m,c_b_n,c_b_nz
real(c_double), intent(in) :: c_b_as
integer(c_size_t), intent(in) :: c_b_ja,c_b_irp
integer(c_size_t), intent(in), value :: c_b_max_row_nz
integer(c_int), intent(in), value :: c_impl_choice
integer(c_size_t), intent(out), value :: c_c_m,c_c_n,c_c_nz
real(c_double), intent(out) :: c_c_as
integer(c_size_t), intent(out) :: c_c_ja,c_c_irp
integer(c_size_t), intent(out), value :: c_c_max_row_nz
integer(psb_ipk_), intent(out) :: c_info
end subroutine psb_f_spmm
end interface
end subroutine dspmm
Loading…
Cancel
Save