xref: /petsc/src/mat/impls/aij/seq/ftn-kernels/fmult.F90 (revision 0113e7197f2fe148981cb35224c0458dc22f510d)
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!
7c96caaccSSatish Balaysubroutine FortranMultTransposeAddAIJ(n,x,ii,jj,a,y)
8c96caaccSSatish Balay  implicit none
9c96caaccSSatish Balay  PetscScalar x(0:*),a(0:*),y(0:*)
10c96caaccSSatish Balay  PetscScalar alpha
11c96caaccSSatish Balay  PetscInt    n,ii(*),jj(0:*)
12c96caaccSSatish Balay
13c96caaccSSatish Balay  PetscInt    i,j,jstart,jend
14c96caaccSSatish Balay
15c96caaccSSatish Balay  jend  = ii(1)
16*0113e719SMartin Diehl  do i=1,n
17c96caaccSSatish Balay    jstart = jend
18c96caaccSSatish Balay    jend   = ii(i+1)
19c96caaccSSatish Balay    alpha  = x(i-1)
20*0113e719SMartin Diehl    do j=jstart,jend-1
21c96caaccSSatish Balay      y(jj(j)) = y(jj(j)) + alpha*a(j)
22*0113e719SMartin Diehl    end do
23*0113e719SMartin Diehl  end do
24c96caaccSSatish Balay
25*0113e719SMartin Diehlend subroutine FortranMultTransposeAddAIJ
26c96caaccSSatish Balay
27c96caaccSSatish Balaysubroutine FortranMultAIJ(n,x,ii,jj,a,y)
28c96caaccSSatish Balay  implicit none
29c96caaccSSatish Balay  PetscScalar x(0:*),a(0:*),y(*)
30c96caaccSSatish Balay  PetscInt    n,ii(*),jj(0:*)
31c96caaccSSatish Balay
32c96caaccSSatish Balay  PetscInt i,j,jstart,jend
33c96caaccSSatish Balay  PetscScalar  sum
34c96caaccSSatish Balay
355161815cSJunchao Zhang#ifdef PETSC_USE_OPENMP_KERNELS
365161815cSJunchao Zhang  !omp parallel do private(j,jstart,jend,sum)
375161815cSJunchao Zhang#endif
38*0113e719SMartin Diehl  do i=1,n
395161815cSJunchao Zhang    jstart = ii(i)
40c96caaccSSatish Balay    jend   = ii(i+1)
41c96caaccSSatish Balay    sum    = 0.d0
42*0113e719SMartin Diehl    do j=jstart,jend-1
43c96caaccSSatish Balay      sum = sum + a(j)*x(jj(j))
44*0113e719SMartin Diehl    end do
45c96caaccSSatish Balay    y(i) = sum
46*0113e719SMartin Diehl  end do
47c96caaccSSatish Balay
48*0113e719SMartin Diehlend subroutine FortranMultAIJ
49