Actual source code: ex79f.F

petsc-3.7.3 2016-07-24
Report Typos and Errors
  1: !
  2: !   This program demonstrates use of MatGetRowIJ() from Fortran
  3: !
  4:       program main
  5:       implicit none
  6: #include <petsc/finclude/petscsys.h>
  7: #include <petsc/finclude/petscmat.h>
  8: #include <petsc/finclude/petscis.h>
  9: #include <petsc/finclude/petscviewer.h>

 11:       Mat         A,Ad,Ao
 12:       PetscErrorCode ierr
 13:       PetscMPIInt rank
 14:       PetscViewer v
 15:       PetscInt i,j,ia(1),ja(1)
 16:       PetscInt n,icol(1),rstart
 17:       PetscInt zero,one,rend
 18:       PetscBool   done
 19:       PetscOffset iia,jja,aaa,iicol
 20:       PetscScalar aa(1)

 22:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 23:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 25:       call PetscViewerBinaryOpen(PETSC_COMM_WORLD,                          &
 26:      & '../../../../share/petsc/datafiles/matrices/' //                       &
 27:      & 'ns-real-int32-float64',                                               &
 28:      &                          FILE_MODE_READ,v,ierr)
 29:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 30:       call MatSetType(A, MATMPIAIJ,ierr)
 31:       call MatLoad(A,v,ierr)
 32:       CHKERRQ(ierr)
 33:       call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr)

 35:       call MatMPIAIJGetSeqAIJ(A,Ad,Ao,icol,iicol,ierr)
 36:       call MatGetOwnershipRange(A,rstart,rend,ierr)
 37: !
 38: !   Print diagonal portion of matrix
 39: !
 40:       call PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1,ierr)
 41:       zero = 0
 42:       one  = 1
 43:       call MatGetRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
 44:       call MatSeqAIJGetArray(Ad,aa,aaa,ierr)
 45:       do 10, i=1,n
 46:         write(7+rank,*) 'row ',i+rstart,' number nonzeros ',                &
 47:      &                   ia(iia+i+1)-ia(iia+i)
 48:         do 20, j=ia(iia+i),ia(iia+i+1)-1
 49:           write(7+rank,*)'  ',j,ja(jja+j)+rstart,aa(aaa+j)
 50:  20     continue
 51:  10   continue
 52:       call MatRestoreRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
 53:       call MatSeqAIJRestoreArray(Ad,aa,aaa,ierr)
 54: !
 55: !   Print off-diagonal portion of matrix
 56: !
 57:       call MatGetRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
 58:       call MatSeqAIJGetArray(Ao,aa,aaa,ierr)
 59:       do 30, i=1,n
 60:         write(7+rank,*) 'row ',i+rstart,' number nonzeros ',               &
 61:      &                  ia(iia+i+1)-ia(iia+i)
 62:         do 40, j=ia(iia+i),ia(iia+i+1)-1
 63:           write(7+rank,*)'  ',j,icol(iicol+ja(jja+j))+1,aa(aaa+j)
 64:  40     continue
 65:  30   continue
 66:       call MatRestoreRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
 67:       call MatSeqAIJRestoreArray(Ao,aa,aaa,ierr)

 69:       call PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1,ierr)

 71:       call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr)
 72:       call MatDestroy(A,ierr)
 73:       call PetscViewerDestroy(v,ierr)
 74:       call PetscFinalize(ierr)
 75:       end