1! 2! 3! Fortran kernel for sparse matrix-vector product in the AIJ matrix format 4! 5#include <petsc/finclude/petscsys.h> 6! 7 subroutine FortranMultTransposeAddAIJ(n,x,ii,jj,a,y) 8 implicit none 9 PetscScalar x(0:*),a(0:*),y(0:*) 10 PetscScalar alpha 11 PetscInt n,ii(*),jj(0:*) 12 13 PetscInt i,j,jstart,jend 14 15 jend = ii(1) 16 do 10,i=1,n 17 jstart = jend 18 jend = ii(i+1) 19 alpha = x(i-1) 20 do 20 j=jstart,jend-1 21 y(jj(j)) = y(jj(j)) + alpha*a(j) 22 20 continue 23 10 continue 24 25 end 26 27 subroutine FortranMultAIJ(n,x,ii,jj,a,y) 28 implicit none 29 PetscScalar x(0:*),a(0:*),y(*) 30 PetscInt n,ii(*),jj(0:*) 31 32 PetscInt i,j,jstart,jend 33 PetscScalar sum 34 35#ifdef PETSC_USE_OPENMP_KERNELS 36 !omp parallel do private(j,jstart,jend,sum) 37#endif 38 do 10,i=1,n 39 jstart = ii(i) 40 jend = ii(i+1) 41 sum = 0.d0 42 do 20 j=jstart,jend-1 43 sum = sum + a(j)*x(jj(j)) 44 20 continue 45 y(i) = sum 46 10 continue 47 48 end 49