xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmultadd.F90 (revision 0ccf82ac36648ce52b79cfc8b55f689a1594b19a)
1c96caaccSSatish Balay!
2c96caaccSSatish Balay!
3c96caaccSSatish Balay!    Fortran kernel for sparse matrix-vector product in the AIJ format
4c96caaccSSatish Balay!
5c96caaccSSatish Balay#include <petsc/finclude/petscsys.h>
6c96caaccSSatish Balay!
7*0ccf82acSMartin Diehlpure subroutine FortranMultAddAIJ(n,x,ii,jj,a,y,z)
8*0ccf82acSMartin Diehl  implicit none (type, external)
9*0ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*),y(*)
10*0ccf82acSMartin Diehl  PetscScalar, intent(inout) :: z(*)
11*0ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
12c96caaccSSatish Balay
13c96caaccSSatish Balay  PetscInt i,j,jstart,jend
14c96caaccSSatish Balay  PetscScalar  sum
15c96caaccSSatish Balay
16c96caaccSSatish Balay  jend  = ii(1)
170113e719SMartin Diehl  do i=1,n
18c96caaccSSatish Balay    jstart = jend
19c96caaccSSatish Balay    jend   = ii(i+1)
20c96caaccSSatish Balay    sum    = y(i)
210113e719SMartin Diehl    do j=jstart,jend-1
22c96caaccSSatish Balay      sum = sum + a(j)*x(jj(j))
230113e719SMartin Diehl    end do
24c96caaccSSatish Balay    z(i) = sum
250113e719SMartin Diehl  end do
260113e719SMartin Diehlend subroutine FortranMultAddAIJ
27