xref: /petsc/src/mat/impls/aij/seq/crl/ftn-kernels/fmultcrl.F90 (revision d66e387e1b17bf1eafe50b0bb7df00ccc9053b5a)
1c96caaccSSatish Balay!
2c96caaccSSatish Balay!
3c96caaccSSatish Balay!    Fortran kernel for sparse matrix-vector product in the AIJ/CRL format
4c96caaccSSatish Balay!
5c96caaccSSatish Balay#include <petsc/finclude/petscsys.h>
6c96caaccSSatish Balay!
70ccf82acSMartin Diehlpure subroutine FortranMultCRL(m,rmax,x,y,icols,acols)
80ccf82acSMartin Diehl  implicit none (type, external)
90ccf82acSMartin Diehl  PetscInt, intent(in) :: m,rmax,icols(m,rmax)
100ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:m-1), acols(m,rmax)
110ccf82acSMartin Diehl  PetscScalar, intent(out) :: y(m)
12c96caaccSSatish Balay
13*d66e387eSMartin Diehl  PetscInt :: i
14c96caaccSSatish Balay
15*d66e387eSMartin Diehl  y(1:m) = acols(1:m,1)*x(icols(1:m,1))
160113e719SMartin Diehl  do i=2,rmax
17*d66e387eSMartin Diehl    y(1:m) = y(1:m) + acols(1:m,i)*x(icols(1:m,i))
180113e719SMartin Diehl  end do
190113e719SMartin Diehlend subroutine FortranMultCRL
20