Actual source code: ex67f.F90
1: !
2: ! This program demonstrates use of MatCreateSubMatrices() from Fortran
3: !
4: program main
5: #include <petsc/finclude/petscmat.h>
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) then
33: nis = 0 ! test nis = 0
34: end if
35: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, nis, zero, PETSC_COPY_VALUES, isrow, ierr))
37: PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_INITIAL_MATRIX, B, ierr))
39: if (rank == 0) then
40: PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
41: end if
43: PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_REUSE_MATRIX, B, ierr))
45: if (rank == 0) then
46: PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
47: end if
49: PetscCallA(ISDestroy(isrow, ierr))
50: PetscCallA(MatDestroy(A, ierr))
51: PetscCallA(MatDestroySubMatrices(nis, B, ierr))
52: PetscCallA(PetscViewerDestroy(v, ierr))
54: PetscCallA(PetscFinalize(ierr))
55: end
57: !/*TEST
58: !
59: ! test:
60: ! requires: double !complex
61: !
62: !TEST*/