xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmultadd.F90 (revision 0ccf82ac36648ce52b79cfc8b55f689a1594b19a)
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  implicit none (type, external)
9  PetscScalar, intent(in) :: x(0:*),a(0:*),y(*)
10  PetscScalar, intent(inout) :: z(*)
11  PetscInt, intent(in) :: n,ii(*),jj(0:*)
12
13  PetscInt i,j,jstart,jend
14  PetscScalar  sum
15
16  jend  = ii(1)
17  do i=1,n
18    jstart = jend
19    jend   = ii(i+1)
20    sum    = y(i)
21    do j=jstart,jend-1
22      sum = sum + a(j)*x(jj(j))
23    end do
24    z(i) = sum
25  end do
26end subroutine FortranMultAddAIJ
27