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