xref: /petsc/src/mat/tests/ex67f.F90 (revision 66af8762ec03dbef0e079729eb2a1734a35ed7ff) !
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
8
9      Mat             A,B(2)
10      PetscErrorCode  ierr
11      PetscInt        nis,zero(1)
12      PetscViewer     v
13      IS              isrow
14      PetscMPIInt     rank
15
16      PetscCallA(PetscInitialize(ierr))
17      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
18
19#if defined(PETSC_USE_64BIT_INDICES)
20      PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,'${PETSC_DIR}/share/petsc/datafiles/matrices/' // 'ns-real-int64-float64',FILE_MODE_READ,v,ierr))
21#else
22      PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,'${PETSC_DIR}/share/petsc/datafiles/matrices/' // 'ns-real-int32-float64',FILE_MODE_READ,v,ierr))
23#endif
24
25      PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
26      PetscCallA(MatSetType(A, MATAIJ,ierr))
27      PetscCallA(MatLoad(A,v,ierr))
28
29      nis     = 1
30      zero(1) = 0
31      if (rank .eq. 1) then
32         nis = 0 ! test nis = 0
33      endif
34      PetscCallA(ISCreateGeneral(PETSC_COMM_SELF,nis,zero,PETSC_COPY_VALUES,isrow,ierr))
35
36      PetscCallA(MatCreateSubmatrices(A,nis,isrow,isrow,MAT_INITIAL_MATRIX,B,ierr))
37
38      if (rank .eq. 0) then
39         PetscCallA(MatView(B(1),PETSC_VIEWER_STDOUT_SELF,ierr))
40      endif
41
42      PetscCallA(MatCreateSubmatrices(A,nis,isrow,isrow,MAT_REUSE_MATRIX,B,ierr))
43
44      if (rank .eq. 0) then
45         PetscCallA(MatView(B(1),PETSC_VIEWER_STDOUT_SELF,ierr))
46      endif
47
48      PetscCallA(ISDestroy(isrow,ierr))
49      PetscCallA(MatDestroy(A,ierr))
50      PetscCallA(MatDestroySubMatrices(nis,B,ierr))
51      PetscCallA(PetscViewerDestroy(v,ierr))
52
53      PetscCallA(PetscFinalize(ierr))
54      end
55
56!/*TEST
57!
58!     test:
59!        requires: double !complex
60!
61!TEST*/
62