xref: /petsc/src/mat/tests/ex63f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1d8606c27SBarry Smith!
2d8606c27SBarry Smith!
3d8606c27SBarry Smith!   This program tests storage of PETSc Dense matrix.
4d8606c27SBarry Smith!   It Creates a Dense matrix, and stores it into a file,
5d8606c27SBarry Smith!   and then reads it back in as a SeqDense and MPIDense
6d8606c27SBarry Smith!   matrix, and prints out the contents.
7d8606c27SBarry Smith!
8d8606c27SBarry Smith#include <petsc/finclude/petscmat.h>
9*c5e229c2SMartin Diehlprogram main
10d8606c27SBarry Smith  use petscmat
11d8606c27SBarry Smith  implicit none
12d8606c27SBarry Smith
13d8606c27SBarry Smith  PetscErrorCode ierr
14d8606c27SBarry Smith  PetscInt row, col, ten
15d8606c27SBarry Smith  PetscMPIInt rank
16d8606c27SBarry Smith  PetscScalar v
17d8606c27SBarry Smith  Mat A
18d8606c27SBarry Smith  PetscViewer view
19d8606c27SBarry Smith
20d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
21d8606c27SBarry Smith
22d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
23d8606c27SBarry Smith!
24d8606c27SBarry Smith!     Proc-0 Create a seq-dense matrix and write it to a file
25d8606c27SBarry Smith!
264820e4eaSBarry Smith  if (rank == 0) then
27d8606c27SBarry Smith    ten = 10
285d83a8b1SBarry Smith    PetscCallA(MatCreateSeqDense(PETSC_COMM_SELF, ten, ten, PETSC_NULL_SCALAR_ARRAY, A, ierr))
29d8606c27SBarry Smith    v = 1.0
30d8606c27SBarry Smith    do row = 0, 9
31d8606c27SBarry Smith      do col = 0, 9
32d8606c27SBarry Smith        PetscCallA(MatSetValue(A, row, col, v, INSERT_VALUES, ierr))
33d8606c27SBarry Smith        v = v + 1.0
34d8606c27SBarry Smith      end do
35d8606c27SBarry Smith    end do
36d8606c27SBarry Smith
37d8606c27SBarry Smith    PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
38d8606c27SBarry Smith    PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
39d8606c27SBarry Smith
40d8606c27SBarry Smith    PetscCallA(PetscObjectSetName(A, 'Original Matrix', ierr))
41d8606c27SBarry Smith    PetscCallA(MatView(A, PETSC_VIEWER_STDOUT_SELF, ierr))
42d8606c27SBarry Smith!
43d8606c27SBarry Smith!        Now Write this matrix to a binary file
44d8606c27SBarry Smith!
45d8606c27SBarry Smith    PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF, 'dense.mat', FILE_MODE_WRITE, view, ierr))
46d8606c27SBarry Smith    PetscCallA(MatView(A, view, ierr))
47d8606c27SBarry Smith    PetscCallA(PetscViewerDestroy(view, ierr))
48d8606c27SBarry Smith    PetscCallA(MatDestroy(A, ierr))
49d8606c27SBarry Smith!
50d8606c27SBarry Smith!        Read this matrix into a SeqDense matrix
51d8606c27SBarry Smith
52d8606c27SBarry Smith    PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF, 'dense.mat', FILE_MODE_READ, view, ierr))
53d8606c27SBarry Smith    PetscCallA(MatCreate(PETSC_COMM_SELF, A, ierr))
54d8606c27SBarry Smith    PetscCallA(MatSetType(A, MATSEQDENSE, ierr))
55d8606c27SBarry Smith    PetscCallA(MatLoad(A, view, ierr))
56d8606c27SBarry Smith
57d8606c27SBarry Smith    PetscCallA(PetscObjectSetName(A, 'SeqDense Matrix read in from file', ierr))
58d8606c27SBarry Smith    PetscCallA(MatView(A, PETSC_VIEWER_STDOUT_SELF, ierr))
59d8606c27SBarry Smith    PetscCallA(MatDestroy(A, ierr))
60d8606c27SBarry Smith    PetscCallA(PetscViewerDestroy(view, ierr))
61d8606c27SBarry Smith  end if
62d8606c27SBarry Smith
63d8606c27SBarry Smith!
64d8606c27SBarry Smith!     Use a barrier, so that the procs do not try opening the file before
65d8606c27SBarry Smith!     it is created.
66d8606c27SBarry Smith!
67d8606c27SBarry Smith  PetscCallMPIA(MPI_Barrier(PETSC_COMM_WORLD, ierr))
68d8606c27SBarry Smith
69d8606c27SBarry Smith  PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD, 'dense.mat', FILE_MODE_READ, view, ierr))
70d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
71d8606c27SBarry Smith  PetscCallA(MatSetType(A, MATMPIDENSE, ierr))
72d8606c27SBarry Smith  PetscCallA(MatLoad(A, view, ierr))
73d8606c27SBarry Smith
74d8606c27SBarry Smith  PetscCallA(PetscObjectSetName(A, 'MPIDense Matrix read in from file', ierr))
75d8606c27SBarry Smith  PetscCallA(MatView(A, PETSC_VIEWER_STDOUT_WORLD, ierr))
76d8606c27SBarry Smith  PetscCallA(MatDestroy(A, ierr))
77d8606c27SBarry Smith  PetscCallA(PetscViewerDestroy(view, ierr))
78d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
79d8606c27SBarry Smithend
80d8606c27SBarry Smith
81d8606c27SBarry Smith!/*TEST
82d8606c27SBarry Smith!
83d8606c27SBarry Smith!   test:
84d8606c27SBarry Smith!      nsize: 2
85d8606c27SBarry Smith!      output_file: output/ex63_1.out
86d8606c27SBarry Smith!
87d8606c27SBarry Smith!TEST*/
88