xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmult.F90 (revision ff45ff59e5033a47c65b211fde43097dd9cb5a3c)
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)
80ccf82acSMartin Diehl  implicit none (type, external)
90ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
100ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(0:*)
110ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
12c96caaccSSatish Balay
13d66e387eSMartin Diehl  PetscInt :: i,jstart,jend
14c96caaccSSatish Balay
15c96caaccSSatish Balay  jend  = ii(1)
160113e719SMartin Diehl  do i=1,n
17c96caaccSSatish Balay    jstart = jend
18c96caaccSSatish Balay    jend   = ii(i+1)
19d66e387eSMartin Diehl    y(jj(jstart:jend-1)) = y(jj(jstart:jend-1)) + x(i-1)*a(jstart:jend-1)
200113e719SMartin Diehl  end do
210113e719SMartin Diehlend subroutine FortranMultTransposeAddAIJ
22c96caaccSSatish Balay
230ccf82acSMartin Diehlpure subroutine FortranMultAIJ(n,x,ii,jj,a,y)
240ccf82acSMartin Diehl  implicit none (type, external)
250ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
260ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(*)
270ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
28c96caaccSSatish Balay
29d66e387eSMartin Diehl  PetscInt :: i,jstart,jend
30c96caaccSSatish Balay
315161815cSJunchao Zhang#ifdef PETSC_USE_OPENMP_KERNELS
32*ff45ff59SMartin Diehl  !omp parallel do private(jstart,jend)
335161815cSJunchao Zhang#endif
340113e719SMartin Diehl  do i=1,n
355161815cSJunchao Zhang    jstart = ii(i)
36c96caaccSSatish Balay    jend   = ii(i+1)
37d66e387eSMartin Diehl    y(i) = sum(a(jstart:jend-1)*x(jj(jstart:jend-1)))
380113e719SMartin Diehl  end do
390113e719SMartin Diehlend subroutine FortranMultAIJ
40