xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmult.F90 (revision fe66ebcc023cb303106674d426ee542bea707d38)
1c96caaccSSatish Balay!
2c96caaccSSatish Balay!
3c96caaccSSatish Balay!    Fortran kernel for sparse matrix-vector product in the AIJ matrix format
4c96caaccSSatish Balay!
5c96caaccSSatish Balay#include <petsc/finclude/petscsys.h>
6c96caaccSSatish Balay!
70ccf82acSMartin Diehlpure subroutine FortranMultTransposeAddAIJ(n,x,ii,jj,a,y)
8*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
90ccf82acSMartin Diehl  implicit none (type, external)
100ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
110ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(0:*)
120ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
13c96caaccSSatish Balay
14d66e387eSMartin Diehl  PetscInt :: i,jstart,jend
15c96caaccSSatish Balay
16c96caaccSSatish Balay  jend = ii(1)
170113e719SMartin Diehl  do i=1,n
18c96caaccSSatish Balay    jstart = jend
19c96caaccSSatish Balay    jend   = ii(i+1)
20d66e387eSMartin Diehl    y(jj(jstart:jend-1)) = y(jj(jstart:jend-1)) + x(i-1)*a(jstart:jend-1)
210113e719SMartin Diehl  end do
220113e719SMartin Diehlend subroutine FortranMultTransposeAddAIJ
23c96caaccSSatish Balay
240ccf82acSMartin Diehlpure subroutine FortranMultAIJ(n,x,ii,jj,a,y)
25*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
260ccf82acSMartin Diehl  implicit none (type, external)
270ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
280ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(*)
290ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
30c96caaccSSatish Balay
31d66e387eSMartin Diehl  PetscInt :: i,jstart,jend
32c96caaccSSatish Balay
335161815cSJunchao Zhang#ifdef PETSC_USE_OPENMP_KERNELS
34ff45ff59SMartin Diehl  !omp parallel do private(jstart,jend)
355161815cSJunchao Zhang#endif
360113e719SMartin Diehl  do i=1,n
375161815cSJunchao Zhang    jstart = ii(i)
38c96caaccSSatish Balay    jend   = ii(i+1)
39d66e387eSMartin Diehl    y(i) = sum(a(jstart:jend-1)*x(jj(jstart:jend-1)))
400113e719SMartin Diehl  end do
410113e719SMartin Diehlend subroutine FortranMultAIJ
42