Actual source code: fmultcrl.F
petsc-3.6.1 2015-08-06
1: !
2: !
3: ! Fortran kernel for sparse matrix-vector product in the AIJ/CRL format
4: !
5: #include <petsc/finclude/petscsysdef.h>
6: !
7: subroutine FortranMultCRL(m,rmax,x,y,icols,acols)
8: implicit none
9: PetscInt m,rmax,icols(m,rmax)
10: PetscScalar x(0:m-1),y(m)
11: PetscScalar acols(m,rmax)
13: PetscInt i,j
15: do 5 j=1,m
16: y(j) = acols(j,1)*x(icols(j,1))
17: 5 continue
19: do 10,i=2,rmax
20: do 20 j=1,m
21: y(j) = y(j) + acols(j,i)*x(icols(j,i))
22: 20 continue
23: 10 continue
25: return
26: end