xref: /petsc/src/mat/tests/ex16f90.F90 (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1!
2!
3!  Tests MatDenseGetArray()
4!
5
6      program main
7#include <petsc/finclude/petscmat.h>
8      use petscmat
9      implicit none
10
11      Mat A
12      PetscErrorCode ierr
13      PetscInt i,j,m,n,iar(1),jar(1)
14      PetscInt one
15      PetscScalar  v(1)
16      PetscScalar, pointer :: array(:,:)
17      PetscMPIInt rank
18      integer :: ashape(2)
19      character(len=80) :: string
20
21      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
22      if (ierr .ne. 0) then
23        print*,'Unable to initialize PETSc'
24        stop
25      endif
26
27      m = 3
28      n = 2
29      one = 1
30!
31!      Create a parallel dense matrix shared by all processors
32!
33      call MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,PETSC_NULL_SCALAR,A,ierr);CHKERRA(ierr)
34
35!
36!     Set values into the matrix. All processors set all values.
37!
38      do 10, i=0,m-1
39        iar(1) = i
40        do 20, j=0,n-1
41          jar(1) = j
42          v(1)   = 9.0/real(i+j+1)
43          call MatSetValues(A,one,iar,one,jar,v,INSERT_VALUES,ierr);CHKERRA(ierr)
44 20     continue
45 10   continue
46
47      call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr);CHKERRA(ierr)
48      call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr);CHKERRA(ierr)
49
50!
51!       Print the matrix to the screen
52!
53      call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr)
54
55!
56!      Print the local matrix shape to the screen for each rank
57!
58      call MatDenseGetArrayF90(A,array,ierr);CHKERRA(ierr)
59      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr);CHKERRA(ierr);
60      ashape = shape(array)
61      write(string, '("[", i0, "]", " shape (", i0, ",", i0, ")", a1)') rank, ashape(1), ashape(2), new_line('a')
62      call PetscSynchronizedPrintf(PETSC_COMM_WORLD, string, ierr);CHKERRA(ierr);
63      call PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT,ierr);CHKERRA(ierr);
64      call MatDenseRestoreArrayF90(A,array,ierr);CHKERRA(ierr)
65!
66!      Free the space used by the matrix
67!
68      call MatDestroy(A,ierr);CHKERRA(ierr)
69      call PetscFinalize(ierr)
70      end
71
72
73!/*TEST
74!
75!   test:
76!      nsize: 2
77!      filter: sort -b
78!      filter_output: sort -b
79!
80!TEST*/
81