Actual source code: fmult.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 matrix format
  4: !
  5: #include <petsc/finclude/petscsysdef.h>
  6: !
  7:       subroutine FortranMultTransposeAddAIJ(n,x,ii,jj,a,y)
  8:       implicit none
  9:       PetscScalar x(0:*),a(0:*),y(0:*)
 10:       PetscScalar alpha
 11:       PetscInt    n,ii(*),jj(0:*)

 13:       PetscInt    i,j,jstart,jend

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

 25:       return
 26:       end

 28:       subroutine FortranMultAIJ(n,x,ii,jj,a,y)
 29:       implicit none
 30:       PetscScalar x(0:*),a(0:*),y(*)
 31:       PetscInt    n,ii(*),jj(0:*)

 33:       PetscInt i,j,jstart,jend
 34:       PetscScalar  sum

 36:       jend  = ii(1)
 37:       do 10,i=1,n
 38:         jstart = jend
 39:         jend   = ii(i+1)
 40:         sum    = 0.d0
 41:         do 20 j=jstart,jend-1
 42:           sum = sum + a(j)*x(jj(j))
 43:  20     continue
 44:         y(i) = sum
 45:  10   continue

 47:       return
 48:       end