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