Actual source code: ex79f.F90

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: !
  2: !   This program demonstrates use of MatGetRowIJ() from Fortran
  3: !
  4:       program main
  5:  #include <petsc/finclude/petscmat.h>
  6:       use petscmat
  7:       implicit none

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

 20:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 21:       if (ierr .ne. 0) then
 22:         print*,'Unable to initialize PETSc'
 23:         stop
 24:       endif
 25:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

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

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

 71:       call PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1,ierr)

 73:       call MatGetDiagonalBlock(A,Ad,ierr)
 74:       call MatView(Ad,PETSC_VIEWER_STDOUT_WORLD,ierr)

 76:       call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr)
 77:       call MatDestroy(A,ierr)
 78:       call PetscViewerDestroy(v,ierr)
 79:       call PetscFinalize(ierr)
 80:       end