xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmultadd.F90 (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1!
2!
3!    Fortran kernel for sparse matrix-vector product in the AIJ format
4!
5#include <petsc/finclude/petscsys.h>
6!
7pure subroutine FortranMultAddAIJ(n,x,ii,jj,a,y,z)
8  use, intrinsic :: ISO_C_binding
9  implicit none (type, external)
10  PetscScalar, intent(in) :: x(0:*),a(0:*),y(*)
11  PetscScalar, intent(inout) :: z(*)
12  PetscInt, intent(in) :: n,ii(*),jj(0:*)
13
14  PetscInt :: i,jstart,jend
15
16  jend  = ii(1)
17  do i=1,n
18    jstart = jend
19    jend   = ii(i+1)
20    z(i) = y(i) + sum(a(jstart:jend-1)*x(jj(jstart:jend-1)))
21  end do
22end subroutine FortranMultAddAIJ
23