Skip to content

Commit

Permalink
sparse linalg: simplify SELLC spmv kernel (#912)
Browse files Browse the repository at this point in the history
  • Loading branch information
perazz authored Jan 3, 2025
2 parents 94f201d + 6590506 commit 05e44f0
Show file tree
Hide file tree
Showing 2 changed files with 13 additions and 13 deletions.
2 changes: 1 addition & 1 deletion src/stdlib_sparse_conversion.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -430,7 +430,7 @@ contains
end block
!-------------------------------------------
! copy values and colum index
allocate(SELLC%col(chunk_size,nnz), source = -1)
allocate(SELLC%col(chunk_size,nnz), source = 1)
allocate(SELLC%data(chunk_size,nnz), source = zero )
block
integer :: lb, ri, iaa, iab, rownnz
Expand Down
24 changes: 12 additions & 12 deletions src/stdlib_sparse_spmv.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -513,23 +513,23 @@ contains
#:for chunk in CHUNKS
pure subroutine chunk_kernel_${chunk}$(n,a,col,x,y)
integer, value :: n
${t1}$, intent(in) :: a(${chunk}$,n), x(*)
${t1}$, intent(in) :: a(${chunk}$,n), x(:)
integer(ilp), intent(in) :: col(${chunk}$,n)
${t1}$, intent(inout) :: y(${chunk}$)
integer :: j
do j = 1, n
where(col(:,j) > 0) y = y + alpha_ * a(:,j) * x(col(:,j))
y(:) = y(:) + alpha_ * a(:,j) * x(col(:,j))
end do
end subroutine
pure subroutine chunk_kernel_trans_${chunk}$(n,a,col,x,y)
integer, value :: n
${t1}$, intent(in) :: a(${chunk}$,n), x(${chunk}$)
integer(ilp), intent(in) :: col(${chunk}$,n)
${t1}$, intent(inout) :: y(*)
${t1}$, intent(inout) :: y(:)
integer :: j, k
do j = 1, n
do k = 1, ${chunk}$
if(col(k,j) > 0) y(col(k,j)) = y(col(k,j)) + alpha_ * a(k,j) * x(k)
y(col(k,j)) = y(col(k,j)) + alpha_ * a(k,j) * x(k)
end do
end do
end subroutine
Expand All @@ -538,11 +538,11 @@ contains
integer, value :: n
${t1}$, intent(in) :: a(${chunk}$,n), x(${chunk}$)
integer(ilp), intent(in) :: col(${chunk}$,n)
${t1}$, intent(inout) :: y(*)
${t1}$, intent(inout) :: y(:)
integer :: j, k
do j = 1, n
do k = 1, ${chunk}$
if(col(k,j) > 0) y(col(k,j)) = y(col(k,j)) + alpha_ * conjg(a(k,j)) * x(k)
y(col(k,j)) = y(col(k,j)) + alpha_ * conjg(a(k,j)) * x(k)
end do
end do
end subroutine
Expand All @@ -551,23 +551,23 @@ contains

pure subroutine chunk_kernel_rm(n,cs,r,a,col,x,y)
integer, value :: n, cs, r
${t1}$, intent(in) :: a(cs,n), x(*)
${t1}$, intent(in) :: a(cs,n), x(:)
integer(ilp), intent(in) :: col(cs,n)
${t1}$, intent(inout) :: y(r)
integer :: j
do j = 1, n
where(col(1:r,j) > 0) y = y + alpha_ * a(1:r,j) * x(col(1:r,j))
y(1:r) = y(1:r) + alpha_ * a(1:r,j) * x(col(1:r,j))
end do
end subroutine
pure subroutine chunk_kernel_rm_trans(n,cs,r,a,col,x,y)
integer, value :: n, cs, r
${t1}$, intent(in) :: a(cs,n), x(r)
integer(ilp), intent(in) :: col(cs,n)
${t1}$, intent(inout) :: y(*)
${t1}$, intent(inout) :: y(:)
integer :: j, k
do j = 1, n
do k = 1, r
if(col(k,j) > 0) y(col(k,j)) = y(col(k,j)) + alpha_ * a(k,j) * x(k)
y(col(k,j)) = y(col(k,j)) + alpha_ * a(k,j) * x(k)
end do
end do
end subroutine
Expand All @@ -576,11 +576,11 @@ contains
integer, value :: n, cs, r
${t1}$, intent(in) :: a(cs,n), x(r)
integer(ilp), intent(in) :: col(cs,n)
${t1}$, intent(inout) :: y(*)
${t1}$, intent(inout) :: y(:)
integer :: j, k
do j = 1, n
do k = 1, r
if(col(k,j) > 0) y(col(k,j)) = y(col(k,j)) + alpha_ * conjg(a(k,j)) * x(k)
y(col(k,j)) = y(col(k,j)) + alpha_ * conjg(a(k,j)) * x(k)
end do
end do
end subroutine
Expand Down

0 comments on commit 05e44f0

Please sign in to comment.