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