Actual source code: fmultadd.F

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: !
  2: !
  3: !    Fortran kernel for sparse matrix-vector product in the AIJ format
  4: !
  5: #include <petsc/finclude/petscsysdef.h>
  6: !
  7:       subroutine FortranMultAddAIJ(n,x,ii,jj,a,y,z)
  8:       implicit none
  9:       PetscScalar      x(0:*),a(0:*),y(*),z(*)
 10:       PetscInt          n,ii(*),jj(0:*)

 12:       PetscInt i,j,jstart,jend
 13:       PetscScalar  sum

 15:       jend  = ii(1)
 16:       do 10,i=1,n
 17:         jstart = jend
 18:         jend   = ii(i+1)
 19:         sum    = y(i)
 20:         do 20 j=jstart,jend-1
 21:           sum = sum + a(j)*x(jj(j))
 22:  20     continue
 23:         z(i) = sum
 24:  10   continue

 26:       return
 27:       end