xref: /petsc/src/mat/tests/ex79f.F90 (revision ce78bad369055609e946c9d2c25ea67a45873e27)
1c4762a1bSJed Brown!
2c4762a1bSJed Brown!   This program demonstrates use of MatGetRowIJ() from Fortran
3c4762a1bSJed Brown!
4c4762a1bSJed Brown      program main
5c4762a1bSJed Brown
6c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
7c4762a1bSJed Brown      use petscmat
8c4762a1bSJed Brown      implicit none
9c4762a1bSJed Brown
10c4762a1bSJed Brown      Mat         A,Ad,Ao
11c4762a1bSJed Brown      PetscErrorCode ierr
12c4762a1bSJed Brown      PetscMPIInt rank
13c4762a1bSJed Brown      PetscViewer v
1442ce371bSBarry Smith      PetscInt i,j
1542ce371bSBarry Smith      PetscInt n,rstart
16c4762a1bSJed Brown      PetscInt zero,one,rend
1742ce371bSBarry Smith      PetscBool   done,bb
1842ce371bSBarry Smith      PetscScalar,pointer :: aa(:)
1942ce371bSBarry Smith      PetscInt,pointer :: ia(:),ja(:),icol(:)
20c4762a1bSJed Brown
21d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
22d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
23c4762a1bSJed Brown
24d8606c27SBarry Smith      PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,'${PETSC_DIR}/share/petsc/datafiles/matrices/' // 'ns-real-int32-float64',FILE_MODE_READ,v,ierr))
25d8606c27SBarry Smith      PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
26d8606c27SBarry Smith      PetscCallA(MatSetType(A, MATMPIAIJ,ierr))
27d8606c27SBarry Smith      PetscCallA(MatLoad(A,v,ierr))
28d8606c27SBarry Smith      PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))
29c4762a1bSJed Brown
30*ce78bad3SBarry Smith      PetscCallA(MatMPIAIJGetSeqAIJ(A,Ad,Ao,icol,ierr))
31d8606c27SBarry Smith      PetscCallA(MatGetOwnershipRange(A,rstart,rend,ierr))
32c4762a1bSJed Brown!
33c4762a1bSJed Brown!   Print diagonal portion of matrix
34c4762a1bSJed Brown!
35d8606c27SBarry Smith      PetscCallA(PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1,ierr))
36c4762a1bSJed Brown      zero = 0
37c4762a1bSJed Brown      one  = 1
3842ce371bSBarry Smith      bb = .true.
39*ce78bad3SBarry Smith      PetscCallA(MatGetRowIJ(Ad,one,bb,bb,n,ia,ja,done,ierr))
40*ce78bad3SBarry Smith      PetscCallA(MatSeqAIJGetArray(Ad,aa,ierr))
41c4762a1bSJed Brown      do 10, i=1,n
4242ce371bSBarry Smith        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',ia(i+1)-ia(i)
4342ce371bSBarry Smith        do 20, j=ia(i),ia(i+1)-1
4442ce371bSBarry Smith          write(7+rank,*)'  ',j,ja(j)+rstart,aa(j)
45c4762a1bSJed Brown 20     continue
46c4762a1bSJed Brown 10   continue
47*ce78bad3SBarry Smith      PetscCallA(MatRestoreRowIJ(Ad,one,bb,bb,n,ia,ja,done,ierr))
48*ce78bad3SBarry Smith      PetscCallA(MatSeqAIJRestoreArray(Ad,aa,ierr))
49c4762a1bSJed Brown!
50c4762a1bSJed Brown!   Print off-diagonal portion of matrix
51c4762a1bSJed Brown!
52*ce78bad3SBarry Smith      PetscCallA(MatGetRowIJ(Ao,one,bb,bb,n,ia,ja,done,ierr))
53*ce78bad3SBarry Smith      PetscCallA(MatSeqAIJGetArray(Ao,aa,ierr))
54c4762a1bSJed Brown      do 30, i=1,n
5542ce371bSBarry Smith        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',ia(i+1)-ia(i)
5642ce371bSBarry Smith        do 40, j=ia(i),ia(i+1)-1
5742ce371bSBarry Smith          write(7+rank,*)'  ',j,icol(ja(j))+1,aa(j)
58c4762a1bSJed Brown 40     continue
59c4762a1bSJed Brown 30   continue
60*ce78bad3SBarry Smith      PetscCallA(MatMPIAIJRestoreSeqAIJ(A,Ad,Ao,icol,ierr))
61*ce78bad3SBarry Smith      PetscCallA(MatRestoreRowIJ(Ao,one,bb,bb,n,ia,ja,done,ierr))
62*ce78bad3SBarry Smith      PetscCallA(MatSeqAIJRestoreArray(Ao,aa,ierr))
63c4762a1bSJed Brown
64d8606c27SBarry Smith      PetscCallA(PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1,ierr))
65c4762a1bSJed Brown
66d8606c27SBarry Smith      PetscCallA(MatGetDiagonalBlock(A,Ad,ierr))
67d8606c27SBarry Smith      PetscCallA(MatView(Ad,PETSC_VIEWER_STDOUT_WORLD,ierr))
68c4762a1bSJed Brown
69d8606c27SBarry Smith      PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))
70d8606c27SBarry Smith      PetscCallA(MatDestroy(A,ierr))
71d8606c27SBarry Smith      PetscCallA(PetscViewerDestroy(v,ierr))
72d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
73c4762a1bSJed Brown      end
74c4762a1bSJed Brown
75c4762a1bSJed Brown!/*TEST
76c4762a1bSJed Brown!
77c4762a1bSJed Brown!     build:
78dfd57a17SPierre Jolivet!       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
79c4762a1bSJed Brown!
80c4762a1bSJed Brown!     test:
81c4762a1bSJed Brown!        args: -binary_read_double -options_left false
82c4762a1bSJed Brown!
83c4762a1bSJed Brown!TEST*/
84