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