xref: /petsc/src/mat/tests/ex79f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
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
14c4762a1bSJed Brown      PetscInt i,j,ia(1),ja(1)
15c4762a1bSJed Brown      PetscInt n,icol(1),rstart
16c4762a1bSJed Brown      PetscInt zero,one,rend
17c4762a1bSJed Brown      PetscBool   done
18c4762a1bSJed Brown      PetscOffset iia,jja,aaa,iicol
19c4762a1bSJed Brown      PetscScalar aa(1)
20c4762a1bSJed Brown
21*d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
22*d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
23c4762a1bSJed Brown
24*d8606c27SBarry Smith      PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,'${PETSC_DIR}/share/petsc/datafiles/matrices/' // 'ns-real-int32-float64',FILE_MODE_READ,v,ierr))
25*d8606c27SBarry Smith      PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
26*d8606c27SBarry Smith      PetscCallA(MatSetType(A, MATMPIAIJ,ierr))
27*d8606c27SBarry Smith      PetscCallA(MatLoad(A,v,ierr))
28*d8606c27SBarry Smith      PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))
29c4762a1bSJed Brown
30*d8606c27SBarry Smith      PetscCallA(MatMPIAIJGetSeqAIJ(A,Ad,Ao,icol,iicol,ierr))
31*d8606c27SBarry Smith      PetscCallA(MatGetOwnershipRange(A,rstart,rend,ierr))
32c4762a1bSJed Brown!
33c4762a1bSJed Brown!   Print diagonal portion of matrix
34c4762a1bSJed Brown!
35*d8606c27SBarry Smith      PetscCallA(PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1,ierr))
36c4762a1bSJed Brown      zero = 0
37c4762a1bSJed Brown      one  = 1
38*d8606c27SBarry Smith      PetscCallA(MatGetRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr))
39*d8606c27SBarry Smith      PetscCallA(MatSeqAIJGetArray(Ad,aa,aaa,ierr))
40c4762a1bSJed Brown      do 10, i=1,n
41*d8606c27SBarry Smith        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',ia(iia+i+1)-ia(iia+i)
42c4762a1bSJed Brown        do 20, j=ia(iia+i),ia(iia+i+1)-1
43c4762a1bSJed Brown          write(7+rank,*)'  ',j,ja(jja+j)+rstart,aa(aaa+j)
44c4762a1bSJed Brown 20     continue
45c4762a1bSJed Brown 10   continue
46*d8606c27SBarry Smith      PetscCallA(MatRestoreRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr))
47*d8606c27SBarry Smith      PetscCallA(MatSeqAIJRestoreArray(Ad,aa,aaa,ierr))
48c4762a1bSJed Brown!
49c4762a1bSJed Brown!   Print off-diagonal portion of matrix
50c4762a1bSJed Brown!
51*d8606c27SBarry Smith      PetscCallA(MatGetRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr))
52*d8606c27SBarry Smith      PetscCallA(MatSeqAIJGetArray(Ao,aa,aaa,ierr))
53c4762a1bSJed Brown      do 30, i=1,n
54*d8606c27SBarry Smith        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',ia(iia+i+1)-ia(iia+i)
55c4762a1bSJed Brown        do 40, j=ia(iia+i),ia(iia+i+1)-1
56c4762a1bSJed Brown          write(7+rank,*)'  ',j,icol(iicol+ja(jja+j))+1,aa(aaa+j)
57c4762a1bSJed Brown 40     continue
58c4762a1bSJed Brown 30   continue
59*d8606c27SBarry Smith      PetscCallA(MatRestoreRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr))
60*d8606c27SBarry Smith      PetscCallA(MatSeqAIJRestoreArray(Ao,aa,aaa,ierr))
61c4762a1bSJed Brown
62*d8606c27SBarry Smith      PetscCallA(PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1,ierr))
63c4762a1bSJed Brown
64*d8606c27SBarry Smith      PetscCallA(MatGetDiagonalBlock(A,Ad,ierr))
65*d8606c27SBarry Smith      PetscCallA(MatView(Ad,PETSC_VIEWER_STDOUT_WORLD,ierr))
66c4762a1bSJed Brown
67*d8606c27SBarry Smith      PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))
68*d8606c27SBarry Smith      PetscCallA(MatDestroy(A,ierr))
69*d8606c27SBarry Smith      PetscCallA(PetscViewerDestroy(v,ierr))
70*d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
71c4762a1bSJed Brown      end
72c4762a1bSJed Brown
73c4762a1bSJed Brown!/*TEST
74c4762a1bSJed Brown!
75c4762a1bSJed Brown!     build:
76dfd57a17SPierre Jolivet!       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
77c4762a1bSJed Brown!
78c4762a1bSJed Brown!     test:
79c4762a1bSJed Brown!        args: -binary_read_double -options_left false
80c4762a1bSJed Brown!
81c4762a1bSJed Brown!TEST*/
82