xref: /petsc/src/mat/tests/ex79f.F90 (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown!
2*c4762a1bSJed Brown!   This program demonstrates use of MatGetRowIJ() from Fortran
3*c4762a1bSJed Brown!
4*c4762a1bSJed Brown      program main
5*c4762a1bSJed Brown
6*c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
7*c4762a1bSJed Brown      use petscmat
8*c4762a1bSJed Brown      implicit none
9*c4762a1bSJed Brown
10*c4762a1bSJed Brown      Mat         A,Ad,Ao
11*c4762a1bSJed Brown      PetscErrorCode ierr
12*c4762a1bSJed Brown      PetscMPIInt rank
13*c4762a1bSJed Brown      PetscViewer v
14*c4762a1bSJed Brown      PetscInt i,j,ia(1),ja(1)
15*c4762a1bSJed Brown      PetscInt n,icol(1),rstart
16*c4762a1bSJed Brown      PetscInt zero,one,rend
17*c4762a1bSJed Brown      PetscBool   done
18*c4762a1bSJed Brown      PetscOffset iia,jja,aaa,iicol
19*c4762a1bSJed Brown      PetscScalar aa(1)
20*c4762a1bSJed Brown
21*c4762a1bSJed Brown      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
22*c4762a1bSJed Brown      if (ierr .ne. 0) then
23*c4762a1bSJed Brown        print*,'Unable to initialize PETSc'
24*c4762a1bSJed Brown        stop
25*c4762a1bSJed Brown      endif
26*c4762a1bSJed Brown      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
27*c4762a1bSJed Brown
28*c4762a1bSJed Brown      call PetscViewerBinaryOpen(PETSC_COMM_WORLD,                          &
29*c4762a1bSJed Brown     & '${PETSC_DIR}/share/petsc/datafiles/matrices/' //                       &
30*c4762a1bSJed Brown     & 'ns-real-int32-float64',                                               &
31*c4762a1bSJed Brown     &                          FILE_MODE_READ,v,ierr)
32*c4762a1bSJed Brown      call MatCreate(PETSC_COMM_WORLD,A,ierr)
33*c4762a1bSJed Brown      call MatSetType(A, MATMPIAIJ,ierr)
34*c4762a1bSJed Brown      call MatLoad(A,v,ierr)
35*c4762a1bSJed Brown      CHKERRA(ierr)
36*c4762a1bSJed Brown      call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr)
37*c4762a1bSJed Brown
38*c4762a1bSJed Brown      call MatMPIAIJGetSeqAIJ(A,Ad,Ao,icol,iicol,ierr)
39*c4762a1bSJed Brown      call MatGetOwnershipRange(A,rstart,rend,ierr)
40*c4762a1bSJed Brown!
41*c4762a1bSJed Brown!   Print diagonal portion of matrix
42*c4762a1bSJed Brown!
43*c4762a1bSJed Brown      call PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1,ierr)
44*c4762a1bSJed Brown      zero = 0
45*c4762a1bSJed Brown      one  = 1
46*c4762a1bSJed Brown      call MatGetRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
47*c4762a1bSJed Brown      call MatSeqAIJGetArray(Ad,aa,aaa,ierr)
48*c4762a1bSJed Brown      do 10, i=1,n
49*c4762a1bSJed Brown        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',                &
50*c4762a1bSJed Brown     &                   ia(iia+i+1)-ia(iia+i)
51*c4762a1bSJed Brown        do 20, j=ia(iia+i),ia(iia+i+1)-1
52*c4762a1bSJed Brown          write(7+rank,*)'  ',j,ja(jja+j)+rstart,aa(aaa+j)
53*c4762a1bSJed Brown 20     continue
54*c4762a1bSJed Brown 10   continue
55*c4762a1bSJed Brown      call MatRestoreRowIJ(Ad,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
56*c4762a1bSJed Brown      call MatSeqAIJRestoreArray(Ad,aa,aaa,ierr)
57*c4762a1bSJed Brown!
58*c4762a1bSJed Brown!   Print off-diagonal portion of matrix
59*c4762a1bSJed Brown!
60*c4762a1bSJed Brown      call MatGetRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
61*c4762a1bSJed Brown      call MatSeqAIJGetArray(Ao,aa,aaa,ierr)
62*c4762a1bSJed Brown      do 30, i=1,n
63*c4762a1bSJed Brown        write(7+rank,*) 'row ',i+rstart,' number nonzeros ',               &
64*c4762a1bSJed Brown     &                  ia(iia+i+1)-ia(iia+i)
65*c4762a1bSJed Brown        do 40, j=ia(iia+i),ia(iia+i+1)-1
66*c4762a1bSJed Brown          write(7+rank,*)'  ',j,icol(iicol+ja(jja+j))+1,aa(aaa+j)
67*c4762a1bSJed Brown 40     continue
68*c4762a1bSJed Brown 30   continue
69*c4762a1bSJed Brown      call MatRestoreRowIJ(Ao,one,zero,zero,n,ia,iia,ja,jja,done,ierr)
70*c4762a1bSJed Brown      call MatSeqAIJRestoreArray(Ao,aa,aaa,ierr)
71*c4762a1bSJed Brown
72*c4762a1bSJed Brown      call PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1,ierr)
73*c4762a1bSJed Brown
74*c4762a1bSJed Brown      call MatGetDiagonalBlock(A,Ad,ierr)
75*c4762a1bSJed Brown      call MatView(Ad,PETSC_VIEWER_STDOUT_WORLD,ierr)
76*c4762a1bSJed Brown
77*c4762a1bSJed Brown      call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr)
78*c4762a1bSJed Brown      call MatDestroy(A,ierr)
79*c4762a1bSJed Brown      call PetscViewerDestroy(v,ierr)
80*c4762a1bSJed Brown      call PetscFinalize(ierr)
81*c4762a1bSJed Brown      end
82*c4762a1bSJed Brown
83*c4762a1bSJed Brown!/*TEST
84*c4762a1bSJed Brown!
85*c4762a1bSJed Brown!     build:
86*c4762a1bSJed Brown!       requires: double !complex !define(PETSC_USE_64BIT_INDICES)
87*c4762a1bSJed Brown!
88*c4762a1bSJed Brown!     test:
89*c4762a1bSJed Brown!        args: -binary_read_double -options_left false
90*c4762a1bSJed Brown!
91*c4762a1bSJed Brown!TEST*/
92