You cannot select more than 25 topics
Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
87 lines
2.3 KiB
Fortran
87 lines
2.3 KiB
Fortran
submodule (psb_d_oacc_hll_mat_mod) psb_d_oacc_hll_csmm_impl
|
|
use psb_base_mod
|
|
contains
|
|
module subroutine psb_d_oacc_hll_csmm(alpha, a, x, beta, y, info, trans)
|
|
implicit none
|
|
class(psb_d_oacc_hll_sparse_mat), intent(in) :: a
|
|
real(psb_dpk_), intent(in) :: alpha, beta
|
|
real(psb_dpk_), intent(in) :: x(:,:)
|
|
real(psb_dpk_), intent(inout) :: y(:,:)
|
|
integer(psb_ipk_), intent(out) :: info
|
|
character, optional, intent(in) :: trans
|
|
|
|
character :: trans_
|
|
integer(psb_ipk_) :: i, j, m, n, k, nxy, nhacks
|
|
logical :: tra
|
|
integer(psb_ipk_) :: err_act
|
|
character(len=20) :: name = 'd_oacc_hll_csmm'
|
|
logical, parameter :: debug = .false.
|
|
|
|
info = psb_success_
|
|
call psb_erractionsave(err_act)
|
|
|
|
if (present(trans)) then
|
|
trans_ = trans
|
|
else
|
|
trans_ = 'N'
|
|
end if
|
|
|
|
if (.not.a%is_asb()) then
|
|
info = psb_err_invalid_mat_state_
|
|
call psb_errpush(info, name)
|
|
goto 9999
|
|
endif
|
|
tra = (psb_toupper(trans_) == 'T') .or. (psb_toupper(trans_) == 'C')
|
|
|
|
if (tra) then
|
|
m = a%get_ncols()
|
|
n = a%get_nrows()
|
|
else
|
|
n = a%get_ncols()
|
|
m = a%get_nrows()
|
|
end if
|
|
|
|
if (size(x,1) < n) then
|
|
info = 36
|
|
call psb_errpush(info, name, i_err = (/3 * ione, n, izero, izero, izero/))
|
|
goto 9999
|
|
end if
|
|
|
|
if (size(y,1) < m) then
|
|
info = 36
|
|
call psb_errpush(info, name, i_err = (/5 * ione, m, izero, izero, izero/))
|
|
goto 9999
|
|
end if
|
|
|
|
if (tra) then
|
|
call a%psb_d_hll_sparse_mat%spmm(alpha, x, beta, y, info, trans)
|
|
else
|
|
nxy = min(size(x,2), size(y,2))
|
|
nhacks = (a%get_nrows() + a%hksz - 1) / a%hksz
|
|
|
|
!$acc parallel loop collapse(2) present(a, x, y)
|
|
do j = 1, nxy
|
|
do i = 1, m
|
|
y(i,j) = beta * y(i,j)
|
|
end do
|
|
end do
|
|
|
|
!$acc parallel loop present(a, x, y)
|
|
do j = 1, nxy
|
|
do k = 1, nhacks
|
|
do i = a%hkoffs(k), a%hkoffs(k + 1) - 1
|
|
y(a%irn(i), j) = y(a%irn(i), j) + alpha * a%val(i) * x(a%ja(i), j)
|
|
end do
|
|
end do
|
|
end do
|
|
endif
|
|
|
|
call psb_erractionrestore(err_act)
|
|
return
|
|
|
|
9999 call psb_error_handler(err_act)
|
|
return
|
|
|
|
end subroutine psb_d_oacc_hll_csmm
|
|
end submodule psb_d_oacc_hll_csmm_impl
|