xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmult.F90 (revision 0ccf82ac36648ce52b79cfc8b55f689a1594b19a)
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!
7*0ccf82acSMartin Diehlpure subroutine FortranMultTransposeAddAIJ(n,x,ii,jj,a,y)
8*0ccf82acSMartin Diehl  implicit none (type, external)
9*0ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
10*0ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(0:*)
11*0ccf82acSMartin Diehl  PetscScalar :: alpha
12*0ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
13c96caaccSSatish Balay
14c96caaccSSatish Balay  PetscInt    i,j,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)
20c96caaccSSatish Balay    alpha  = x(i-1)
210113e719SMartin Diehl    do j=jstart,jend-1
22c96caaccSSatish Balay      y(jj(j)) = y(jj(j)) + alpha*a(j)
230113e719SMartin Diehl    end do
240113e719SMartin Diehl  end do
250113e719SMartin Diehlend subroutine FortranMultTransposeAddAIJ
26c96caaccSSatish Balay
27*0ccf82acSMartin Diehlpure subroutine FortranMultAIJ(n,x,ii,jj,a,y)
28*0ccf82acSMartin Diehl  implicit none (type, external)
29*0ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(0:*),a(0:*)
30*0ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(*)
31*0ccf82acSMartin Diehl  PetscInt, intent(in) :: n,ii(*),jj(0:*)
32c96caaccSSatish Balay
33c96caaccSSatish Balay  PetscInt i,j,jstart,jend
34c96caaccSSatish Balay  PetscScalar  sum
35c96caaccSSatish Balay
365161815cSJunchao Zhang#ifdef PETSC_USE_OPENMP_KERNELS
375161815cSJunchao Zhang  !omp parallel do private(j,jstart,jend,sum)
385161815cSJunchao Zhang#endif
390113e719SMartin Diehl  do i=1,n
405161815cSJunchao Zhang    jstart = ii(i)
41c96caaccSSatish Balay    jend   = ii(i+1)
42c96caaccSSatish Balay    sum    = 0.d0
430113e719SMartin Diehl    do j=jstart,jend-1
44c96caaccSSatish Balay      sum = sum + a(j)*x(jj(j))
450113e719SMartin Diehl    end do
46c96caaccSSatish Balay    y(i) = sum
470113e719SMartin Diehl  end do
480113e719SMartin Diehlend subroutine FortranMultAIJ
49