xref: /petsc/src/mat/tests/ex67f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1d8606c27SBarry Smith!
2d8606c27SBarry Smith!   This program demonstrates use of MatCreateSubMatrices() from Fortran
3d8606c27SBarry Smith!
4d8606c27SBarry Smith#include <petsc/finclude/petscmat.h>
5*c5e229c2SMartin Diehlprogram main
6d8606c27SBarry Smith  use petscmat
7d8606c27SBarry Smith  implicit none
8d8606c27SBarry Smith
9ce78bad3SBarry Smith  Mat A
10ce78bad3SBarry Smith  Mat, pointer :: B(:)
11d8606c27SBarry Smith  PetscErrorCode ierr
12d8606c27SBarry Smith  PetscInt nis, zero(1)
13d8606c27SBarry Smith  PetscViewer v
14d8606c27SBarry Smith  IS isrow
15d8606c27SBarry Smith  PetscMPIInt rank
16d8606c27SBarry Smith
17d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
18d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
19d8606c27SBarry Smith
20d8606c27SBarry Smith#if defined(PETSC_USE_64BIT_INDICES)
21d8606c27SBarry Smith  PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD, '${PETSC_DIR}/share/petsc/datafiles/matrices/'//'ns-real-int64-float64', FILE_MODE_READ, v, ierr))
22d8606c27SBarry Smith#else
23d8606c27SBarry Smith  PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD, '${PETSC_DIR}/share/petsc/datafiles/matrices/'//'ns-real-int32-float64', FILE_MODE_READ, v, ierr))
24d8606c27SBarry Smith#endif
25d8606c27SBarry Smith
26d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
27d8606c27SBarry Smith  PetscCallA(MatSetType(A, MATAIJ, ierr))
28d8606c27SBarry Smith  PetscCallA(MatLoad(A, v, ierr))
29d8606c27SBarry Smith
30d8606c27SBarry Smith  nis = 1
31d8606c27SBarry Smith  zero(1) = 0
324820e4eaSBarry Smith  if (rank == 1) then
33d8606c27SBarry Smith    nis = 0 ! test nis = 0
34d8606c27SBarry Smith  end if
35d8606c27SBarry Smith  PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, nis, zero, PETSC_COPY_VALUES, isrow, ierr))
36d8606c27SBarry Smith
37ce78bad3SBarry Smith  PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_INITIAL_MATRIX, B, ierr))
38d8606c27SBarry Smith
394820e4eaSBarry Smith  if (rank == 0) then
40d8606c27SBarry Smith    PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
41d8606c27SBarry Smith  end if
42d8606c27SBarry Smith
43ce78bad3SBarry Smith  PetscCallA(MatCreateSubmatrices(A, nis, [isrow], [isrow], MAT_REUSE_MATRIX, B, ierr))
44d8606c27SBarry Smith
454820e4eaSBarry Smith  if (rank == 0) then
46d8606c27SBarry Smith    PetscCallA(MatView(B(1), PETSC_VIEWER_STDOUT_SELF, ierr))
47d8606c27SBarry Smith  end if
48d8606c27SBarry Smith
49d8606c27SBarry Smith  PetscCallA(ISDestroy(isrow, ierr))
50d8606c27SBarry Smith  PetscCallA(MatDestroy(A, ierr))
51d8606c27SBarry Smith  PetscCallA(MatDestroySubMatrices(nis, B, ierr))
52d8606c27SBarry Smith  PetscCallA(PetscViewerDestroy(v, ierr))
53d8606c27SBarry Smith
54d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
55d8606c27SBarry Smithend
56d8606c27SBarry Smith
57d8606c27SBarry Smith!/*TEST
58d8606c27SBarry Smith!
59d8606c27SBarry Smith!     test:
60d8606c27SBarry Smith!        requires: double !complex
61d8606c27SBarry Smith!
62d8606c27SBarry Smith!TEST*/
63