Actual source code: ex67f.F90

  1: !
  2: !   This program demonstrates use of MatCreateSubMatrices() from Fortran
  3: !
  4: #include <petsc/finclude/petscmat.h>
  5: program main
  6:   use petscmat
  7:   implicit none

  9:   Mat A
 10:   Mat, pointer :: B(:)
 11:   PetscErrorCode ierr
 12:   PetscInt nis, zero(1)
 13:   PetscViewer v
 14:   IS isrow
 15:   PetscMPIInt rank

 17:   PetscCallA(PetscInitialize(ierr))
 18:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))

 20: #if defined(PETSC_USE_64BIT_INDICES)
 21:   PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD, '${PETSC_DIR}/share/petsc/datafiles/matrices/'//'ns-real-int64-float64', FILE_MODE_READ, v, ierr))
 22: #else
 23:   PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD, '${PETSC_DIR}/share/petsc/datafiles/matrices/'//'ns-real-int32-float64', FILE_MODE_READ, v, ierr))
 24: #endif

 26:   PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
 27:   PetscCallA(MatSetType(A, MATAIJ, ierr))
 28:   PetscCallA(MatLoad(A, v, ierr))

 30:   nis = 1
 31:   zero(1) = 0
 32:   if (rank == 1) nis = 0 ! test nis = 0
 33:   PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, nis, zero, PETSC_COPY_VALUES, isrow, ierr))

 35:   PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_INITIAL_MATRIX, B, ierr))

 37:   if (rank == 0) then
 38:     PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
 39:   end if

 41:   PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_REUSE_MATRIX, B, ierr))

 43:   if (rank == 0) then
 44:     PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
 45:   end if

 47:   PetscCallA(ISDestroy(isrow, ierr))
 48:   PetscCallA(MatDestroy(A, ierr))
 49:   PetscCallA(MatDestroySubMatrices(nis, B, ierr))
 50:   PetscCallA(PetscViewerDestroy(v, ierr))

 52:   PetscCallA(PetscFinalize(ierr))
 53: end

 55: !/*TEST
 56: !
 57: !     test:
 58: !        requires: double !complex
 59: !
 60: !TEST*/